xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 36bce6e8a8185c8c15b8a5a56b7d3ed7fe3f62ed)
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);
26a297a907SKarl Rupp   if (hdf5->file_id) 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   herr_t            herr;
745c6c1daeSBarry Smith   PetscErrorCode    ierr;
755c6c1daeSBarry Smith 
765c6c1daeSBarry Smith   PetscFunctionBegin;
775c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
785c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
795c6c1daeSBarry Smith   plist_id = H5Pcreate(H5P_FILE_ACCESS);
805c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
81ce94432eSBarry Smith   herr = H5Pset_fapl_mpio(plist_id, PetscObjectComm((PetscObject)viewer), info);CHKERRQ(herr);
825c6c1daeSBarry Smith #endif
835c6c1daeSBarry Smith   /* Create or open the file collectively */
845c6c1daeSBarry Smith   switch (hdf5->btype) {
855c6c1daeSBarry Smith   case FILE_MODE_READ:
865c6c1daeSBarry Smith     hdf5->file_id = H5Fopen(name, H5F_ACC_RDONLY, plist_id);
875c6c1daeSBarry Smith     break;
885c6c1daeSBarry Smith   case FILE_MODE_APPEND:
895c6c1daeSBarry Smith     hdf5->file_id = H5Fopen(name, H5F_ACC_RDWR, plist_id);
905c6c1daeSBarry Smith     break;
915c6c1daeSBarry Smith   case FILE_MODE_WRITE:
925c6c1daeSBarry Smith     hdf5->file_id = H5Fcreate(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id);
935c6c1daeSBarry Smith     break;
945c6c1daeSBarry Smith   default:
955c6c1daeSBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
965c6c1daeSBarry Smith   }
975c6c1daeSBarry Smith   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
985c6c1daeSBarry Smith   H5Pclose(plist_id);
995c6c1daeSBarry Smith   PetscFunctionReturn(0);
1005c6c1daeSBarry Smith }
1015c6c1daeSBarry Smith 
1025c6c1daeSBarry Smith #undef __FUNCT__
1035c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerCreate_HDF5"
1048cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
1055c6c1daeSBarry Smith {
1065c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
1075c6c1daeSBarry Smith   PetscErrorCode   ierr;
1085c6c1daeSBarry Smith 
1095c6c1daeSBarry Smith   PetscFunctionBegin;
110b00a9115SJed Brown   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
1115c6c1daeSBarry Smith 
1125c6c1daeSBarry Smith   v->data         = (void*) hdf5;
1135c6c1daeSBarry Smith   v->ops->destroy = PetscViewerDestroy_HDF5;
1145c6c1daeSBarry Smith   v->ops->flush   = 0;
1155c6c1daeSBarry Smith   hdf5->btype     = (PetscFileMode) -1;
1165c6c1daeSBarry Smith   hdf5->filename  = 0;
1175c6c1daeSBarry Smith   hdf5->timestep  = -1;
1180298fd71SBarry Smith   hdf5->groups    = NULL;
1195c6c1daeSBarry Smith 
1200b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
1210b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
1225c6c1daeSBarry Smith   PetscFunctionReturn(0);
1235c6c1daeSBarry Smith }
1245c6c1daeSBarry Smith 
1255c6c1daeSBarry Smith #undef __FUNCT__
1265c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5Open"
1275c6c1daeSBarry Smith /*@C
1285c6c1daeSBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
1295c6c1daeSBarry Smith 
1305c6c1daeSBarry Smith    Collective on MPI_Comm
1315c6c1daeSBarry Smith 
1325c6c1daeSBarry Smith    Input Parameters:
1335c6c1daeSBarry Smith +  comm - MPI communicator
1345c6c1daeSBarry Smith .  name - name of file
1355c6c1daeSBarry Smith -  type - type of file
1365c6c1daeSBarry Smith $    FILE_MODE_WRITE - create new file for binary output
1375c6c1daeSBarry Smith $    FILE_MODE_READ - open existing file for binary input
1385c6c1daeSBarry Smith $    FILE_MODE_APPEND - open existing file for binary output
1395c6c1daeSBarry Smith 
1405c6c1daeSBarry Smith    Output Parameter:
1415c6c1daeSBarry Smith .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
1425c6c1daeSBarry Smith 
1435c6c1daeSBarry Smith    Level: beginner
1445c6c1daeSBarry Smith 
1455c6c1daeSBarry Smith    Note:
1465c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
1475c6c1daeSBarry Smith 
1485c6c1daeSBarry Smith    Concepts: HDF5 files
1495c6c1daeSBarry Smith    Concepts: PetscViewerHDF5^creating
1505c6c1daeSBarry Smith 
1515c6c1daeSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerSetFormat(), PetscViewerDestroy(),
1525c6c1daeSBarry Smith           VecView(), MatView(), VecLoad(), MatLoad(),
1535c6c1daeSBarry Smith           PetscFileMode, PetscViewer
1545c6c1daeSBarry Smith @*/
1555c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
1565c6c1daeSBarry Smith {
1575c6c1daeSBarry Smith   PetscErrorCode ierr;
1585c6c1daeSBarry Smith 
1595c6c1daeSBarry Smith   PetscFunctionBegin;
1605c6c1daeSBarry Smith   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
1615c6c1daeSBarry Smith   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
1625c6c1daeSBarry Smith   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
1635c6c1daeSBarry Smith   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
1645c6c1daeSBarry Smith   PetscFunctionReturn(0);
1655c6c1daeSBarry Smith }
1665c6c1daeSBarry Smith 
1675c6c1daeSBarry Smith #undef __FUNCT__
1685c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5GetFileId"
1695c6c1daeSBarry Smith /*@C
1705c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
1715c6c1daeSBarry Smith 
1725c6c1daeSBarry Smith   Not collective
1735c6c1daeSBarry Smith 
1745c6c1daeSBarry Smith   Input Parameter:
1755c6c1daeSBarry Smith . viewer - the PetscViewer
1765c6c1daeSBarry Smith 
1775c6c1daeSBarry Smith   Output Parameter:
1785c6c1daeSBarry Smith . file_id - The file id
1795c6c1daeSBarry Smith 
1805c6c1daeSBarry Smith   Level: intermediate
1815c6c1daeSBarry Smith 
1825c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open()
1835c6c1daeSBarry Smith @*/
1845c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
1855c6c1daeSBarry Smith {
1865c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1875c6c1daeSBarry Smith 
1885c6c1daeSBarry Smith   PetscFunctionBegin;
1895c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1905c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
1915c6c1daeSBarry Smith   PetscFunctionReturn(0);
1925c6c1daeSBarry Smith }
1935c6c1daeSBarry Smith 
1945c6c1daeSBarry Smith #undef __FUNCT__
1955c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5PushGroup"
1965c6c1daeSBarry Smith /*@C
1975c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
1985c6c1daeSBarry Smith 
1995c6c1daeSBarry Smith   Not collective
2005c6c1daeSBarry Smith 
2015c6c1daeSBarry Smith   Input Parameters:
2025c6c1daeSBarry Smith + viewer - the PetscViewer
2035c6c1daeSBarry Smith - name - The group name
2045c6c1daeSBarry Smith 
2055c6c1daeSBarry Smith   Level: intermediate
2065c6c1daeSBarry Smith 
2075c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
2085c6c1daeSBarry Smith @*/
2095c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
2105c6c1daeSBarry Smith {
2115c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2125c6c1daeSBarry Smith   GroupList        *groupNode;
2135c6c1daeSBarry Smith   PetscErrorCode   ierr;
2145c6c1daeSBarry Smith 
2155c6c1daeSBarry Smith   PetscFunctionBegin;
2165c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2175c6c1daeSBarry Smith   PetscValidCharPointer(name,2);
2185c6c1daeSBarry Smith   ierr = PetscMalloc(sizeof(GroupList), &groupNode);CHKERRQ(ierr);
2195c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
220a297a907SKarl Rupp 
2215c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
2225c6c1daeSBarry Smith   hdf5->groups    = groupNode;
2235c6c1daeSBarry Smith   PetscFunctionReturn(0);
2245c6c1daeSBarry Smith }
2255c6c1daeSBarry Smith 
2265c6c1daeSBarry Smith #undef __FUNCT__
2275c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5PopGroup"
2283ef9c667SSatish Balay /*@
2295c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
2305c6c1daeSBarry Smith 
2315c6c1daeSBarry Smith   Not collective
2325c6c1daeSBarry Smith 
2335c6c1daeSBarry Smith   Input Parameter:
2345c6c1daeSBarry Smith . viewer - the PetscViewer
2355c6c1daeSBarry Smith 
2365c6c1daeSBarry Smith   Level: intermediate
2375c6c1daeSBarry Smith 
2385c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup()
2395c6c1daeSBarry Smith @*/
2405c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
2415c6c1daeSBarry Smith {
2425c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2435c6c1daeSBarry Smith   GroupList        *groupNode;
2445c6c1daeSBarry Smith   PetscErrorCode   ierr;
2455c6c1daeSBarry Smith 
2465c6c1daeSBarry Smith   PetscFunctionBegin;
2475c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
24882f516ccSBarry Smith   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
2495c6c1daeSBarry Smith   groupNode    = hdf5->groups;
2505c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
2515c6c1daeSBarry Smith   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
2525c6c1daeSBarry Smith   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
2535c6c1daeSBarry Smith   PetscFunctionReturn(0);
2545c6c1daeSBarry Smith }
2555c6c1daeSBarry Smith 
2565c6c1daeSBarry Smith #undef __FUNCT__
2575c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5GetGroup"
2585c6c1daeSBarry Smith /*@C
2590298fd71SBarry Smith   PetscViewerHDF5GetGroup - Get the current HDF5 group for output. If none has been assigned, returns NULL.
2605c6c1daeSBarry Smith 
2615c6c1daeSBarry Smith   Not collective
2625c6c1daeSBarry Smith 
2635c6c1daeSBarry Smith   Input Parameter:
2645c6c1daeSBarry Smith . viewer - the PetscViewer
2655c6c1daeSBarry Smith 
2665c6c1daeSBarry Smith   Output Parameter:
2675c6c1daeSBarry Smith . name - The group name
2685c6c1daeSBarry Smith 
2695c6c1daeSBarry Smith   Level: intermediate
2705c6c1daeSBarry Smith 
2715c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup()
2725c6c1daeSBarry Smith @*/
2735c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
2745c6c1daeSBarry Smith {
2755c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
2765c6c1daeSBarry Smith 
2775c6c1daeSBarry Smith   PetscFunctionBegin;
2785c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2795c6c1daeSBarry Smith   PetscValidCharPointer(name,2);
280a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
2810298fd71SBarry Smith   else *name = NULL;
2825c6c1daeSBarry Smith   PetscFunctionReturn(0);
2835c6c1daeSBarry Smith }
2845c6c1daeSBarry Smith 
2855c6c1daeSBarry Smith #undef __FUNCT__
2865c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5IncrementTimestep"
2873ef9c667SSatish Balay /*@
2885c6c1daeSBarry Smith   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
2895c6c1daeSBarry Smith 
2905c6c1daeSBarry Smith   Not collective
2915c6c1daeSBarry Smith 
2925c6c1daeSBarry Smith   Input Parameter:
2935c6c1daeSBarry Smith . viewer - the PetscViewer
2945c6c1daeSBarry Smith 
2955c6c1daeSBarry Smith   Level: intermediate
2965c6c1daeSBarry Smith 
2975c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
2985c6c1daeSBarry Smith @*/
2995c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
3005c6c1daeSBarry Smith {
3015c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3025c6c1daeSBarry Smith 
3035c6c1daeSBarry Smith   PetscFunctionBegin;
3045c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
3055c6c1daeSBarry Smith   ++hdf5->timestep;
3065c6c1daeSBarry Smith   PetscFunctionReturn(0);
3075c6c1daeSBarry Smith }
3085c6c1daeSBarry Smith 
3095c6c1daeSBarry Smith #undef __FUNCT__
3105c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5SetTimestep"
3113ef9c667SSatish Balay /*@
3125c6c1daeSBarry Smith   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
3135c6c1daeSBarry Smith   of -1 disables blocking with timesteps.
3145c6c1daeSBarry Smith 
3155c6c1daeSBarry Smith   Not collective
3165c6c1daeSBarry Smith 
3175c6c1daeSBarry Smith   Input Parameters:
3185c6c1daeSBarry Smith + viewer - the PetscViewer
3195c6c1daeSBarry Smith - timestep - The timestep number
3205c6c1daeSBarry Smith 
3215c6c1daeSBarry Smith   Level: intermediate
3225c6c1daeSBarry Smith 
3235c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
3245c6c1daeSBarry Smith @*/
3255c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
3265c6c1daeSBarry Smith {
3275c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3285c6c1daeSBarry Smith 
3295c6c1daeSBarry Smith   PetscFunctionBegin;
3305c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
3315c6c1daeSBarry Smith   hdf5->timestep = timestep;
3325c6c1daeSBarry Smith   PetscFunctionReturn(0);
3335c6c1daeSBarry Smith }
3345c6c1daeSBarry Smith 
3355c6c1daeSBarry Smith #undef __FUNCT__
3365c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5GetTimestep"
3373ef9c667SSatish Balay /*@
3385c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
3395c6c1daeSBarry Smith 
3405c6c1daeSBarry Smith   Not collective
3415c6c1daeSBarry Smith 
3425c6c1daeSBarry Smith   Input Parameter:
3435c6c1daeSBarry Smith . viewer - the PetscViewer
3445c6c1daeSBarry Smith 
3455c6c1daeSBarry Smith   Output Parameter:
3465c6c1daeSBarry Smith . timestep - The timestep number
3475c6c1daeSBarry Smith 
3485c6c1daeSBarry Smith   Level: intermediate
3495c6c1daeSBarry Smith 
3505c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
3515c6c1daeSBarry Smith @*/
3525c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
3535c6c1daeSBarry Smith {
3545c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3555c6c1daeSBarry Smith 
3565c6c1daeSBarry Smith   PetscFunctionBegin;
3575c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
3585c6c1daeSBarry Smith   PetscValidPointer(timestep,2);
3595c6c1daeSBarry Smith   *timestep = hdf5->timestep;
3605c6c1daeSBarry Smith   PetscFunctionReturn(0);
3615c6c1daeSBarry Smith }
3625c6c1daeSBarry Smith 
363a75e6a4aSMatthew G. Knepley 
364*36bce6e8SMatthew G. Knepley #undef __FUNCT__
365*36bce6e8SMatthew G. Knepley #define __FUNCT__ "PetscDataTypeToHDF5DataType"
366*36bce6e8SMatthew G. Knepley /*@C
367*36bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
368*36bce6e8SMatthew G. Knepley 
369*36bce6e8SMatthew G. Knepley   Not collective
370*36bce6e8SMatthew G. Knepley 
371*36bce6e8SMatthew G. Knepley   Input Parameter:
372*36bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
373*36bce6e8SMatthew G. Knepley 
374*36bce6e8SMatthew G. Knepley   Output Parameter:
375*36bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
376*36bce6e8SMatthew G. Knepley 
377*36bce6e8SMatthew G. Knepley   Level: advanced
378*36bce6e8SMatthew G. Knepley 
379*36bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
380*36bce6e8SMatthew G. Knepley @*/
381*36bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
382*36bce6e8SMatthew G. Knepley {
383*36bce6e8SMatthew G. Knepley   PetscFunctionBegin;
384*36bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
385*36bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
386*36bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_LLONG;
387*36bce6e8SMatthew G. Knepley #else
388*36bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_INT;
389*36bce6e8SMatthew G. Knepley #endif
390*36bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
391*36bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
392*36bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
393*36bce6e8SMatthew G. Knepley   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_DOUBLE;
394*36bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_DOUBLE;
395*36bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
396*36bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
397*36bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
398*36bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
399*36bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
400*36bce6e8SMatthew G. Knepley }
401*36bce6e8SMatthew G. Knepley 
402*36bce6e8SMatthew G. Knepley #undef __FUNCT__
403*36bce6e8SMatthew G. Knepley #define __FUNCT__ "PetscHDF5DataTypeToPetscDataType"
404*36bce6e8SMatthew G. Knepley /*@C
405*36bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
406*36bce6e8SMatthew G. Knepley 
407*36bce6e8SMatthew G. Knepley   Not collective
408*36bce6e8SMatthew G. Knepley 
409*36bce6e8SMatthew G. Knepley   Input Parameter:
410*36bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
411*36bce6e8SMatthew G. Knepley 
412*36bce6e8SMatthew G. Knepley   Output Parameter:
413*36bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
414*36bce6e8SMatthew G. Knepley 
415*36bce6e8SMatthew G. Knepley   Level: advanced
416*36bce6e8SMatthew G. Knepley 
417*36bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
418*36bce6e8SMatthew G. Knepley @*/
419*36bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
420*36bce6e8SMatthew G. Knepley {
421*36bce6e8SMatthew G. Knepley   PetscFunctionBegin;
422*36bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
423*36bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
424*36bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
425*36bce6e8SMatthew G. Knepley #else
426*36bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
427*36bce6e8SMatthew G. Knepley #endif
428*36bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
429*36bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
430*36bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
431*36bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
432*36bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
433*36bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
434*36bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
435*36bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
436*36bce6e8SMatthew G. Knepley }
437*36bce6e8SMatthew G. Knepley 
438*36bce6e8SMatthew G. Knepley #undef __FUNCT__
439*36bce6e8SMatthew G. Knepley #define __FUNCT__ "PetscViewerHDF5WriteAttribute"
440*36bce6e8SMatthew G. Knepley /*@
441*36bce6e8SMatthew G. Knepley  PetscViewerHDF5WriteAttribute - Write a scalar attribute
442*36bce6e8SMatthew G. Knepley 
443*36bce6e8SMatthew G. Knepley   Input Parameters:
444*36bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer
445*36bce6e8SMatthew G. Knepley . parent - The parent name
446*36bce6e8SMatthew G. Knepley . name   - The attribute name
447*36bce6e8SMatthew G. Knepley . datatype - The attribute type
448*36bce6e8SMatthew G. Knepley - value    - The attribute value
449*36bce6e8SMatthew G. Knepley 
450*36bce6e8SMatthew G. Knepley   Level: advanced
451*36bce6e8SMatthew G. Knepley 
452*36bce6e8SMatthew G. Knepley .seealso: PetscViewerHDF5Open()
453*36bce6e8SMatthew G. Knepley @*/
454*36bce6e8SMatthew G. Knepley PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
455*36bce6e8SMatthew G. Knepley {
456*36bce6e8SMatthew G. Knepley   hid_t          h5, dataspace, dataset, attribute, dtype, status;
457*36bce6e8SMatthew G. Knepley   PetscErrorCode ierr;
458*36bce6e8SMatthew G. Knepley 
459*36bce6e8SMatthew G. Knepley   PetscFunctionBegin;
460*36bce6e8SMatthew G. Knepley   PetscValidPointer(parent, 2);
461*36bce6e8SMatthew G. Knepley   PetscValidPointer(name, 3);
462*36bce6e8SMatthew G. Knepley   PetscValidPointer(value, 4);
463*36bce6e8SMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
464*36bce6e8SMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
465*36bce6e8SMatthew G. Knepley   dataspace = H5Screate(H5S_SCALAR);
466*36bce6e8SMatthew G. Knepley   if (dataspace < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_LIB, "Could not create dataspace for attribute %s of %s", name, parent);
467*36bce6e8SMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
468*36bce6e8SMatthew G. Knepley   dataset = H5Dopen2(h5, parent, H5P_DEFAULT);
469*36bce6e8SMatthew G. Knepley   if (dataset < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_LIB, "Could not open parent dataset for attribute %s of %s", name, parent);
470*36bce6e8SMatthew G. Knepley   attribute = H5Acreate2(dataset, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT);
471*36bce6e8SMatthew G. Knepley #else
472*36bce6e8SMatthew G. Knepley   dataset = H5Dopen(h5, parent);
473*36bce6e8SMatthew G. Knepley   if (dataset < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_LIB, "Could not open parent dataset for attribute %s of %s", name, parent);
474*36bce6e8SMatthew G. Knepley   attribute = H5Acreate(dataset, name, dtype, dataspace, H5P_DEFAULT);
475*36bce6e8SMatthew G. Knepley #endif
476*36bce6e8SMatthew G. Knepley   if (attribute < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_LIB, "Could not create attribute %s of %s", name, parent);
477*36bce6e8SMatthew G. Knepley   status = H5Awrite(attribute, dtype, value);CHKERRQ(status);
478*36bce6e8SMatthew G. Knepley   status = H5Aclose(attribute);CHKERRQ(status);
479*36bce6e8SMatthew G. Knepley   status = H5Dclose(dataset);CHKERRQ(status);
480*36bce6e8SMatthew G. Knepley   status = H5Sclose(dataspace);CHKERRQ(status);
481*36bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
482*36bce6e8SMatthew G. Knepley }
483*36bce6e8SMatthew G. Knepley 
484a75e6a4aSMatthew G. Knepley /*
485a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
486a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
487a75e6a4aSMatthew G. Knepley */
488a75e6a4aSMatthew G. Knepley static int Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
489a75e6a4aSMatthew G. Knepley 
490a75e6a4aSMatthew G. Knepley #undef __FUNCT__
491a75e6a4aSMatthew G. Knepley #define __FUNCT__ "PETSC_VIEWER_HDF5_"
492a75e6a4aSMatthew G. Knepley /*@C
493a75e6a4aSMatthew G. Knepley   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
494a75e6a4aSMatthew G. Knepley 
495a75e6a4aSMatthew G. Knepley   Collective on MPI_Comm
496a75e6a4aSMatthew G. Knepley 
497a75e6a4aSMatthew G. Knepley   Input Parameter:
498a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer
499a75e6a4aSMatthew G. Knepley 
500a75e6a4aSMatthew G. Knepley   Level: intermediate
501a75e6a4aSMatthew G. Knepley 
502a75e6a4aSMatthew G. Knepley   Options Database Keys:
503a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name>
504a75e6a4aSMatthew G. Knepley 
505a75e6a4aSMatthew G. Knepley   Environmental variables:
506a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME
507a75e6a4aSMatthew G. Knepley 
508a75e6a4aSMatthew G. Knepley   Notes:
509a75e6a4aSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
510a75e6a4aSMatthew G. Knepley   an error code.  The HDF5 PetscViewer is usually used in the form
511a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
512a75e6a4aSMatthew G. Knepley 
513a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
514a75e6a4aSMatthew G. Knepley @*/
515a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
516a75e6a4aSMatthew G. Knepley {
517a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
518a75e6a4aSMatthew G. Knepley   PetscBool      flg;
519a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
520a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
521a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
522a75e6a4aSMatthew G. Knepley 
523a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
524a75e6a4aSMatthew 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);}
525a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
526a75e6a4aSMatthew G. Knepley     ierr = MPI_Keyval_create(MPI_NULL_COPY_FN,MPI_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
527a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
528a75e6a4aSMatthew G. Knepley   }
529a75e6a4aSMatthew G. Knepley   ierr = MPI_Attr_get(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
530a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
531a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
532a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
533a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
534a75e6a4aSMatthew G. Knepley     if (!flg) {
535a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname,"output.h5");
536a75e6a4aSMatthew G. Knepley       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
537a75e6a4aSMatthew G. Knepley     }
538a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
539a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
540a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
541a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
542a75e6a4aSMatthew G. Knepley     ierr = MPI_Attr_put(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
543a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
544a75e6a4aSMatthew G. Knepley   }
545a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
546a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
547a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
548a75e6a4aSMatthew G. Knepley }
549*36bce6e8SMatthew G. Knepley 
5505c6c1daeSBarry Smith #if defined(oldhdf4stuff)
5515c6c1daeSBarry Smith #undef __FUNCT__
5525c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5WriteSDS"
5535c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5WriteSDS(PetscViewer viewer, float *xf, int d, int *dims,int bs)
5545c6c1daeSBarry Smith {
5555c6c1daeSBarry Smith   int              i;
5565c6c1daeSBarry Smith   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
5575c6c1daeSBarry Smith   int32            sds_id,zero32[3],dims32[3];
5585c6c1daeSBarry Smith 
5595c6c1daeSBarry Smith   PetscFunctionBegin;
5605c6c1daeSBarry Smith   for (i = 0; i < d; i++) {
5615c6c1daeSBarry Smith     zero32[i] = 0;
5625c6c1daeSBarry Smith     dims32[i] = dims[i];
5635c6c1daeSBarry Smith   }
5645c6c1daeSBarry Smith   sds_id = SDcreate(vhdf5->sd_id, "Vec", DFNT_FLOAT32, d, dims32);
5655c6c1daeSBarry Smith   if (sds_id < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SDcreate failed");
5665c6c1daeSBarry Smith   SDwritedata(sds_id, zero32, 0, dims32, xf);
5675c6c1daeSBarry Smith   SDendaccess(sds_id);
5685c6c1daeSBarry Smith   PetscFunctionReturn(0);
5695c6c1daeSBarry Smith }
5705c6c1daeSBarry Smith 
5715c6c1daeSBarry Smith #endif
572