xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 6056859292bf77bcb84dfd2567a78b2e9798b702)
1af0996ceSBarry 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;
1582ea9b62SBarry Smith   PetscBool     basedimension2;  /* save vectors and DMDA vectors with a dimension of at least 2 even if the bs/dof is 1 */
169a0502c6SHåkon Strandenes   PetscBool     spoutput;  /* write data in single precision even if PETSc is compiled with double precision PetscReal */
175c6c1daeSBarry Smith } PetscViewer_HDF5;
185c6c1daeSBarry Smith 
194416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v)
2082ea9b62SBarry Smith {
2182ea9b62SBarry Smith   PetscErrorCode   ierr;
2282ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;
2382ea9b62SBarry Smith 
2482ea9b62SBarry Smith   PetscFunctionBegin;
2582ea9b62SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr);
2682ea9b62SBarry Smith   ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr);
279a0502c6SHåkon Strandenes   ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr);
2882ea9b62SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
2982ea9b62SBarry Smith   PetscFunctionReturn(0);
3082ea9b62SBarry Smith }
3182ea9b62SBarry Smith 
325c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_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 = PetscFree(hdf5->filename);CHKERRQ(ierr);
39729ab6d9SBarry Smith   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
405c6c1daeSBarry Smith   PetscFunctionReturn(0);
415c6c1daeSBarry Smith }
425c6c1daeSBarry Smith 
435c6c1daeSBarry Smith PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
445c6c1daeSBarry Smith {
455c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
465c6c1daeSBarry Smith   PetscErrorCode   ierr;
475c6c1daeSBarry Smith 
485c6c1daeSBarry Smith   PetscFunctionBegin;
495c6c1daeSBarry Smith   ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr);
505c6c1daeSBarry Smith   while (hdf5->groups) {
515c6c1daeSBarry Smith     GroupList *tmp = hdf5->groups->next;
525c6c1daeSBarry Smith 
535c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups->name);CHKERRQ(ierr);
545c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups);CHKERRQ(ierr);
555c6c1daeSBarry Smith     hdf5->groups = tmp;
565c6c1daeSBarry Smith   }
575c6c1daeSBarry Smith   ierr = PetscFree(hdf5);CHKERRQ(ierr);
580b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
59d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr);
600b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
615c6c1daeSBarry Smith   PetscFunctionReturn(0);
625c6c1daeSBarry Smith }
635c6c1daeSBarry Smith 
645c6c1daeSBarry Smith PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
655c6c1daeSBarry Smith {
665c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
675c6c1daeSBarry Smith 
685c6c1daeSBarry Smith   PetscFunctionBegin;
695c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
705c6c1daeSBarry Smith   hdf5->btype = type;
715c6c1daeSBarry Smith   PetscFunctionReturn(0);
725c6c1daeSBarry Smith }
735c6c1daeSBarry Smith 
7482ea9b62SBarry Smith PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
7582ea9b62SBarry Smith {
7682ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7782ea9b62SBarry Smith 
7882ea9b62SBarry Smith   PetscFunctionBegin;
7982ea9b62SBarry Smith   hdf5->basedimension2 = flg;
8082ea9b62SBarry Smith   PetscFunctionReturn(0);
8182ea9b62SBarry Smith }
8282ea9b62SBarry Smith 
83df863907SAlex Fikl /*@
8482ea9b62SBarry Smith      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
8582ea9b62SBarry Smith        dimension of 2.
8682ea9b62SBarry Smith 
8782ea9b62SBarry Smith     Logically Collective on PetscViewer
8882ea9b62SBarry Smith 
8982ea9b62SBarry Smith   Input Parameters:
9082ea9b62SBarry Smith +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
9182ea9b62SBarry Smith -  flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1
9282ea9b62SBarry Smith 
9382ea9b62SBarry Smith   Options Database:
9482ea9b62SBarry Smith .  -viewer_hdf5_base_dimension2 - turns on (true) or off (false) using a dimension of 2 in the HDF5 file even if the bs/dof of the vector is 1
9582ea9b62SBarry Smith 
9682ea9b62SBarry Smith 
9782ea9b62SBarry Smith   Notes: Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
9882ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
9982ea9b62SBarry Smith 
10082ea9b62SBarry Smith   Level: intermediate
10182ea9b62SBarry Smith 
10282ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
10382ea9b62SBarry Smith 
10482ea9b62SBarry Smith @*/
10582ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
10682ea9b62SBarry Smith {
10782ea9b62SBarry Smith   PetscErrorCode ierr;
10882ea9b62SBarry Smith 
10982ea9b62SBarry Smith   PetscFunctionBegin;
11082ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
11182ea9b62SBarry Smith   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
11282ea9b62SBarry Smith   PetscFunctionReturn(0);
11382ea9b62SBarry Smith }
11482ea9b62SBarry Smith 
115df863907SAlex Fikl /*@
11682ea9b62SBarry Smith      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
11782ea9b62SBarry Smith        dimension of 2.
11882ea9b62SBarry Smith 
11982ea9b62SBarry Smith     Logically Collective on PetscViewer
12082ea9b62SBarry Smith 
12182ea9b62SBarry Smith   Input Parameter:
12282ea9b62SBarry Smith .  viewer - the PetscViewer, must be of type HDF5
12382ea9b62SBarry Smith 
12482ea9b62SBarry Smith   Output Parameter:
12582ea9b62SBarry Smith .  flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1
12682ea9b62SBarry Smith 
12782ea9b62SBarry Smith   Notes: Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
12882ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
12982ea9b62SBarry Smith 
13082ea9b62SBarry Smith   Level: intermediate
13182ea9b62SBarry Smith 
13282ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
13382ea9b62SBarry Smith 
13482ea9b62SBarry Smith @*/
13582ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
13682ea9b62SBarry Smith {
13782ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
13882ea9b62SBarry Smith 
13982ea9b62SBarry Smith   PetscFunctionBegin;
14082ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
14182ea9b62SBarry Smith   *flg = hdf5->basedimension2;
14282ea9b62SBarry Smith   PetscFunctionReturn(0);
14382ea9b62SBarry Smith }
14482ea9b62SBarry Smith 
1459a0502c6SHåkon Strandenes PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
1469a0502c6SHåkon Strandenes {
1479a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1489a0502c6SHåkon Strandenes 
1499a0502c6SHåkon Strandenes   PetscFunctionBegin;
1509a0502c6SHåkon Strandenes   hdf5->spoutput = flg;
1519a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
1529a0502c6SHåkon Strandenes }
1539a0502c6SHåkon Strandenes 
154df863907SAlex Fikl /*@
1559a0502c6SHåkon Strandenes      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
1569a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
1579a0502c6SHåkon Strandenes 
1589a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
1599a0502c6SHåkon Strandenes 
1609a0502c6SHåkon Strandenes   Input Parameters:
1619a0502c6SHåkon Strandenes +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
1629a0502c6SHåkon Strandenes -  flg - if PETSC_TRUE the data will be written to disk with single precision
1639a0502c6SHåkon Strandenes 
1649a0502c6SHåkon Strandenes   Options Database:
1659a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
1669a0502c6SHåkon Strandenes 
1679a0502c6SHåkon Strandenes 
1689a0502c6SHåkon Strandenes   Notes: Setting this option does not make any difference if PETSc is compiled with single precision
1699a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
1709a0502c6SHåkon Strandenes 
1719a0502c6SHåkon Strandenes   Level: intermediate
1729a0502c6SHåkon Strandenes 
1739a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
1749a0502c6SHåkon Strandenes           PetscReal
1759a0502c6SHåkon Strandenes 
1769a0502c6SHåkon Strandenes @*/
1779a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
1789a0502c6SHåkon Strandenes {
1799a0502c6SHåkon Strandenes   PetscErrorCode ierr;
1809a0502c6SHåkon Strandenes 
1819a0502c6SHåkon Strandenes   PetscFunctionBegin;
1829a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1839a0502c6SHåkon Strandenes   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
1849a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
1859a0502c6SHåkon Strandenes }
1869a0502c6SHåkon Strandenes 
187df863907SAlex Fikl /*@
1889a0502c6SHåkon Strandenes      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
1899a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
1909a0502c6SHåkon Strandenes 
1919a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
1929a0502c6SHåkon Strandenes 
1939a0502c6SHåkon Strandenes   Input Parameter:
1949a0502c6SHåkon Strandenes .  viewer - the PetscViewer, must be of type HDF5
1959a0502c6SHåkon Strandenes 
1969a0502c6SHåkon Strandenes   Output Parameter:
1979a0502c6SHåkon Strandenes .  flg - if PETSC_TRUE the data will be written to disk with single precision
1989a0502c6SHåkon Strandenes 
1999a0502c6SHåkon Strandenes   Notes: Setting this option does not make any difference if PETSc is compiled with single precision
2009a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2019a0502c6SHåkon Strandenes 
2029a0502c6SHåkon Strandenes   Level: intermediate
2039a0502c6SHåkon Strandenes 
2049a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2059a0502c6SHåkon Strandenes           PetscReal
2069a0502c6SHåkon Strandenes 
2079a0502c6SHåkon Strandenes @*/
2089a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
2099a0502c6SHåkon Strandenes {
2109a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2119a0502c6SHåkon Strandenes 
2129a0502c6SHåkon Strandenes   PetscFunctionBegin;
2139a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2149a0502c6SHåkon Strandenes   *flg = hdf5->spoutput;
2159a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2169a0502c6SHåkon Strandenes }
2179a0502c6SHåkon Strandenes 
2185c6c1daeSBarry Smith PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
2195c6c1daeSBarry Smith {
2205c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2215c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
2225c6c1daeSBarry Smith   MPI_Info          info = MPI_INFO_NULL;
2235c6c1daeSBarry Smith #endif
2245c6c1daeSBarry Smith   hid_t             plist_id;
2255c6c1daeSBarry Smith   PetscErrorCode    ierr;
2265c6c1daeSBarry Smith 
2275c6c1daeSBarry Smith   PetscFunctionBegin;
228f683cc58SToby Isaac   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
229f683cc58SToby Isaac   if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);}
2305c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
2315c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
232729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
2335c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
234729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info));
2355c6c1daeSBarry Smith #endif
2365c6c1daeSBarry Smith   /* Create or open the file collectively */
2375c6c1daeSBarry Smith   switch (hdf5->btype) {
2385c6c1daeSBarry Smith   case FILE_MODE_READ:
239729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
2405c6c1daeSBarry Smith     break;
2415c6c1daeSBarry Smith   case FILE_MODE_APPEND:
242729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
2435c6c1daeSBarry Smith     break;
2445c6c1daeSBarry Smith   case FILE_MODE_WRITE:
245729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
2465c6c1daeSBarry Smith     break;
2475c6c1daeSBarry Smith   default:
2485c6c1daeSBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
2495c6c1daeSBarry Smith   }
2505c6c1daeSBarry Smith   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
251729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
2525c6c1daeSBarry Smith   PetscFunctionReturn(0);
2535c6c1daeSBarry Smith }
2545c6c1daeSBarry Smith 
255d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
256d1232d7fSToby Isaac {
257d1232d7fSToby Isaac   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
258d1232d7fSToby Isaac 
259d1232d7fSToby Isaac   PetscFunctionBegin;
260d1232d7fSToby Isaac   *name = vhdf5->filename;
261d1232d7fSToby Isaac   PetscFunctionReturn(0);
262d1232d7fSToby Isaac }
263d1232d7fSToby Isaac 
2648556b5ebSBarry Smith /*MC
2658556b5ebSBarry Smith    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
2668556b5ebSBarry Smith 
2678556b5ebSBarry Smith 
2688556b5ebSBarry Smith .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
2698556b5ebSBarry Smith            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
2708556b5ebSBarry Smith            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
2718556b5ebSBarry Smith            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
2728556b5ebSBarry Smith 
2731b266c99SBarry Smith   Level: beginner
2748556b5ebSBarry Smith M*/
275d1232d7fSToby Isaac 
2768cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
2775c6c1daeSBarry Smith {
2785c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
2795c6c1daeSBarry Smith   PetscErrorCode   ierr;
2805c6c1daeSBarry Smith 
2815c6c1daeSBarry Smith   PetscFunctionBegin;
282b00a9115SJed Brown   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
2835c6c1daeSBarry Smith 
2845c6c1daeSBarry Smith   v->data                = (void*) hdf5;
2855c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
28682ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
2875c6c1daeSBarry Smith   v->ops->flush          = 0;
2885c6c1daeSBarry Smith   hdf5->btype            = (PetscFileMode) -1;
2895c6c1daeSBarry Smith   hdf5->filename         = 0;
2905c6c1daeSBarry Smith   hdf5->timestep         = -1;
2910298fd71SBarry Smith   hdf5->groups           = NULL;
2925c6c1daeSBarry Smith 
2930b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
294d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr);
2950b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
29682ea9b62SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
2979a0502c6SHåkon Strandenes   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr);
2985c6c1daeSBarry Smith   PetscFunctionReturn(0);
2995c6c1daeSBarry Smith }
3005c6c1daeSBarry Smith 
3015c6c1daeSBarry Smith /*@C
3025c6c1daeSBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
3035c6c1daeSBarry Smith 
3045c6c1daeSBarry Smith    Collective on MPI_Comm
3055c6c1daeSBarry Smith 
3065c6c1daeSBarry Smith    Input Parameters:
3075c6c1daeSBarry Smith +  comm - MPI communicator
3085c6c1daeSBarry Smith .  name - name of file
3095c6c1daeSBarry Smith -  type - type of file
3105c6c1daeSBarry Smith $    FILE_MODE_WRITE - create new file for binary output
3115c6c1daeSBarry Smith $    FILE_MODE_READ - open existing file for binary input
3125c6c1daeSBarry Smith $    FILE_MODE_APPEND - open existing file for binary output
3135c6c1daeSBarry Smith 
3145c6c1daeSBarry Smith    Output Parameter:
3155c6c1daeSBarry Smith .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
3165c6c1daeSBarry Smith 
31782ea9b62SBarry Smith   Options Database:
31882ea9b62SBarry Smith .  -viewer_hdf5_base_dimension2 - turns on (true) or off (false) using a dimension of 2 in the HDF5 file even if the bs/dof of the vector is 1
3199a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
32082ea9b62SBarry Smith 
3215c6c1daeSBarry Smith    Level: beginner
3225c6c1daeSBarry Smith 
3235c6c1daeSBarry Smith    Note:
3245c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
3255c6c1daeSBarry Smith 
3265c6c1daeSBarry Smith    Concepts: HDF5 files
3275c6c1daeSBarry Smith    Concepts: PetscViewerHDF5^creating
3285c6c1daeSBarry Smith 
3296a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
3309a0502c6SHåkon Strandenes           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
3319a0502c6SHåkon Strandenes           MatLoad(), PetscFileMode, PetscViewer
3325c6c1daeSBarry Smith @*/
3335c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
3345c6c1daeSBarry Smith {
3355c6c1daeSBarry Smith   PetscErrorCode ierr;
3365c6c1daeSBarry Smith 
3375c6c1daeSBarry Smith   PetscFunctionBegin;
3385c6c1daeSBarry Smith   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
3395c6c1daeSBarry Smith   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
3405c6c1daeSBarry Smith   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
3415c6c1daeSBarry Smith   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
3425c6c1daeSBarry Smith   PetscFunctionReturn(0);
3435c6c1daeSBarry Smith }
3445c6c1daeSBarry Smith 
3455c6c1daeSBarry Smith /*@C
3465c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
3475c6c1daeSBarry Smith 
3485c6c1daeSBarry Smith   Not collective
3495c6c1daeSBarry Smith 
3505c6c1daeSBarry Smith   Input Parameter:
3515c6c1daeSBarry Smith . viewer - the PetscViewer
3525c6c1daeSBarry Smith 
3535c6c1daeSBarry Smith   Output Parameter:
3545c6c1daeSBarry Smith . file_id - The file id
3555c6c1daeSBarry Smith 
3565c6c1daeSBarry Smith   Level: intermediate
3575c6c1daeSBarry Smith 
3585c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open()
3595c6c1daeSBarry Smith @*/
3605c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
3615c6c1daeSBarry Smith {
3625c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3635c6c1daeSBarry Smith 
3645c6c1daeSBarry Smith   PetscFunctionBegin;
3655c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
3665c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
3675c6c1daeSBarry Smith   PetscFunctionReturn(0);
3685c6c1daeSBarry Smith }
3695c6c1daeSBarry Smith 
3705c6c1daeSBarry Smith /*@C
3715c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
3725c6c1daeSBarry Smith 
3735c6c1daeSBarry Smith   Not collective
3745c6c1daeSBarry Smith 
3755c6c1daeSBarry Smith   Input Parameters:
3765c6c1daeSBarry Smith + viewer - the PetscViewer
3775c6c1daeSBarry Smith - name - The group name
3785c6c1daeSBarry Smith 
3795c6c1daeSBarry Smith   Level: intermediate
3805c6c1daeSBarry Smith 
3815c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
3825c6c1daeSBarry Smith @*/
3835c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
3845c6c1daeSBarry Smith {
3855c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3865c6c1daeSBarry Smith   GroupList        *groupNode;
3875c6c1daeSBarry Smith   PetscErrorCode   ierr;
3885c6c1daeSBarry Smith 
3895c6c1daeSBarry Smith   PetscFunctionBegin;
3905c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
3915c6c1daeSBarry Smith   PetscValidCharPointer(name,2);
39295dccacaSBarry Smith   ierr = PetscNew(&groupNode);CHKERRQ(ierr);
3935c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
394a297a907SKarl Rupp 
3955c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
3965c6c1daeSBarry Smith   hdf5->groups    = groupNode;
3975c6c1daeSBarry Smith   PetscFunctionReturn(0);
3985c6c1daeSBarry Smith }
3995c6c1daeSBarry Smith 
4003ef9c667SSatish Balay /*@
4015c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
4025c6c1daeSBarry Smith 
4035c6c1daeSBarry Smith   Not collective
4045c6c1daeSBarry Smith 
4055c6c1daeSBarry Smith   Input Parameter:
4065c6c1daeSBarry Smith . viewer - the PetscViewer
4075c6c1daeSBarry Smith 
4085c6c1daeSBarry Smith   Level: intermediate
4095c6c1daeSBarry Smith 
4105c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup()
4115c6c1daeSBarry Smith @*/
4125c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
4135c6c1daeSBarry Smith {
4145c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4155c6c1daeSBarry Smith   GroupList        *groupNode;
4165c6c1daeSBarry Smith   PetscErrorCode   ierr;
4175c6c1daeSBarry Smith 
4185c6c1daeSBarry Smith   PetscFunctionBegin;
4195c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
42082f516ccSBarry Smith   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
4215c6c1daeSBarry Smith   groupNode    = hdf5->groups;
4225c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
4235c6c1daeSBarry Smith   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
4245c6c1daeSBarry Smith   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
4255c6c1daeSBarry Smith   PetscFunctionReturn(0);
4265c6c1daeSBarry Smith }
4275c6c1daeSBarry Smith 
4285c6c1daeSBarry Smith /*@C
4290298fd71SBarry Smith   PetscViewerHDF5GetGroup - Get the current HDF5 group for output. If none has been assigned, returns NULL.
4305c6c1daeSBarry Smith 
4315c6c1daeSBarry Smith   Not collective
4325c6c1daeSBarry Smith 
4335c6c1daeSBarry Smith   Input Parameter:
4345c6c1daeSBarry Smith . viewer - the PetscViewer
4355c6c1daeSBarry Smith 
4365c6c1daeSBarry Smith   Output Parameter:
4375c6c1daeSBarry Smith . name - The group name
4385c6c1daeSBarry Smith 
4395c6c1daeSBarry Smith   Level: intermediate
4405c6c1daeSBarry Smith 
4415c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup()
4425c6c1daeSBarry Smith @*/
4435c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
4445c6c1daeSBarry Smith {
4455c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
4465c6c1daeSBarry Smith 
4475c6c1daeSBarry Smith   PetscFunctionBegin;
4485c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
449c959eef4SJed Brown   PetscValidPointer(name,2);
450a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
4510298fd71SBarry Smith   else *name = NULL;
4525c6c1daeSBarry Smith   PetscFunctionReturn(0);
4535c6c1daeSBarry Smith }
4545c6c1daeSBarry Smith 
4553ef9c667SSatish Balay /*@
4565c6c1daeSBarry Smith   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
4575c6c1daeSBarry Smith 
4585c6c1daeSBarry Smith   Not collective
4595c6c1daeSBarry Smith 
4605c6c1daeSBarry Smith   Input Parameter:
4615c6c1daeSBarry Smith . viewer - the PetscViewer
4625c6c1daeSBarry Smith 
4635c6c1daeSBarry Smith   Level: intermediate
4645c6c1daeSBarry Smith 
4655c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
4665c6c1daeSBarry Smith @*/
4675c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
4685c6c1daeSBarry Smith {
4695c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4705c6c1daeSBarry Smith 
4715c6c1daeSBarry Smith   PetscFunctionBegin;
4725c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
4735c6c1daeSBarry Smith   ++hdf5->timestep;
4745c6c1daeSBarry Smith   PetscFunctionReturn(0);
4755c6c1daeSBarry Smith }
4765c6c1daeSBarry Smith 
4773ef9c667SSatish Balay /*@
4785c6c1daeSBarry Smith   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
4795c6c1daeSBarry Smith   of -1 disables blocking with timesteps.
4805c6c1daeSBarry Smith 
4815c6c1daeSBarry Smith   Not collective
4825c6c1daeSBarry Smith 
4835c6c1daeSBarry Smith   Input Parameters:
4845c6c1daeSBarry Smith + viewer - the PetscViewer
4855c6c1daeSBarry Smith - timestep - The timestep number
4865c6c1daeSBarry Smith 
4875c6c1daeSBarry Smith   Level: intermediate
4885c6c1daeSBarry Smith 
4895c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
4905c6c1daeSBarry Smith @*/
4915c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
4925c6c1daeSBarry Smith {
4935c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4945c6c1daeSBarry Smith 
4955c6c1daeSBarry Smith   PetscFunctionBegin;
4965c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
4975c6c1daeSBarry Smith   hdf5->timestep = timestep;
4985c6c1daeSBarry Smith   PetscFunctionReturn(0);
4995c6c1daeSBarry Smith }
5005c6c1daeSBarry Smith 
5013ef9c667SSatish Balay /*@
5025c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
5035c6c1daeSBarry Smith 
5045c6c1daeSBarry Smith   Not collective
5055c6c1daeSBarry Smith 
5065c6c1daeSBarry Smith   Input Parameter:
5075c6c1daeSBarry Smith . viewer - the PetscViewer
5085c6c1daeSBarry Smith 
5095c6c1daeSBarry Smith   Output Parameter:
5105c6c1daeSBarry Smith . timestep - The timestep number
5115c6c1daeSBarry Smith 
5125c6c1daeSBarry Smith   Level: intermediate
5135c6c1daeSBarry Smith 
5145c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
5155c6c1daeSBarry Smith @*/
5165c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
5175c6c1daeSBarry Smith {
5185c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5195c6c1daeSBarry Smith 
5205c6c1daeSBarry Smith   PetscFunctionBegin;
5215c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5225c6c1daeSBarry Smith   PetscValidPointer(timestep,2);
5235c6c1daeSBarry Smith   *timestep = hdf5->timestep;
5245c6c1daeSBarry Smith   PetscFunctionReturn(0);
5255c6c1daeSBarry Smith }
5265c6c1daeSBarry Smith 
52736bce6e8SMatthew G. Knepley /*@C
52836bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
52936bce6e8SMatthew G. Knepley 
53036bce6e8SMatthew G. Knepley   Not collective
53136bce6e8SMatthew G. Knepley 
53236bce6e8SMatthew G. Knepley   Input Parameter:
53336bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
53436bce6e8SMatthew G. Knepley 
53536bce6e8SMatthew G. Knepley   Output Parameter:
53636bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
53736bce6e8SMatthew G. Knepley 
53836bce6e8SMatthew G. Knepley   Level: advanced
53936bce6e8SMatthew G. Knepley 
54036bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
54136bce6e8SMatthew G. Knepley @*/
54236bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
54336bce6e8SMatthew G. Knepley {
54436bce6e8SMatthew G. Knepley   PetscFunctionBegin;
54536bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
54636bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
54736bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_LLONG;
54836bce6e8SMatthew G. Knepley #else
54936bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_INT;
55036bce6e8SMatthew G. Knepley #endif
55136bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
55236bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
55336bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
55436bce6e8SMatthew G. Knepley   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_DOUBLE;
55536bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_DOUBLE;
55636bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
55736bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
55836bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
5597e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
56036bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
56136bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
56236bce6e8SMatthew G. Knepley }
56336bce6e8SMatthew G. Knepley 
56436bce6e8SMatthew G. Knepley /*@C
56536bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
56636bce6e8SMatthew G. Knepley 
56736bce6e8SMatthew G. Knepley   Not collective
56836bce6e8SMatthew G. Knepley 
56936bce6e8SMatthew G. Knepley   Input Parameter:
57036bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
57136bce6e8SMatthew G. Knepley 
57236bce6e8SMatthew G. Knepley   Output Parameter:
57336bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
57436bce6e8SMatthew G. Knepley 
57536bce6e8SMatthew G. Knepley   Level: advanced
57636bce6e8SMatthew G. Knepley 
57736bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
57836bce6e8SMatthew G. Knepley @*/
57936bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
58036bce6e8SMatthew G. Knepley {
58136bce6e8SMatthew G. Knepley   PetscFunctionBegin;
58236bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
58336bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
58436bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
58536bce6e8SMatthew G. Knepley #else
58636bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
58736bce6e8SMatthew G. Knepley #endif
58836bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
58936bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
59036bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
59136bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
59236bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
59336bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
5947e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
59536bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
59636bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
59736bce6e8SMatthew G. Knepley }
59836bce6e8SMatthew G. Knepley 
599df863907SAlex Fikl /*@C
60036bce6e8SMatthew G. Knepley  PetscViewerHDF5WriteAttribute - Write a scalar attribute
60136bce6e8SMatthew G. Knepley 
60236bce6e8SMatthew G. Knepley   Input Parameters:
60336bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer
60436bce6e8SMatthew G. Knepley . parent - The parent name
60536bce6e8SMatthew G. Knepley . name   - The attribute name
60636bce6e8SMatthew G. Knepley . datatype - The attribute type
60736bce6e8SMatthew G. Knepley - value    - The attribute value
60836bce6e8SMatthew G. Knepley 
60936bce6e8SMatthew G. Knepley   Level: advanced
61036bce6e8SMatthew G. Knepley 
611e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute()
61236bce6e8SMatthew G. Knepley @*/
61336bce6e8SMatthew G. Knepley PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
61436bce6e8SMatthew G. Knepley {
615*60568592SMatthew G. Knepley   hid_t          h5, dataspace, obj, attribute, dtype;
61636bce6e8SMatthew G. Knepley   PetscErrorCode ierr;
61736bce6e8SMatthew G. Knepley 
61836bce6e8SMatthew G. Knepley   PetscFunctionBegin;
6195cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
62036bce6e8SMatthew G. Knepley   PetscValidPointer(parent, 2);
62136bce6e8SMatthew G. Knepley   PetscValidPointer(name, 3);
62236bce6e8SMatthew G. Knepley   PetscValidPointer(value, 4);
62336bce6e8SMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
6247e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
6257e97c476SMatthew G. Knepley     size_t len;
6267e97c476SMatthew G. Knepley     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
627729ab6d9SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
6287e97c476SMatthew G. Knepley   }
62936bce6e8SMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
630729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
631*60568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
63236bce6e8SMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
633*60568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
63436bce6e8SMatthew G. Knepley #else
635*60568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Acreate,(obj, name, dtype, dataspace, H5P_DEFAULT));
63636bce6e8SMatthew G. Knepley #endif
637729ab6d9SBarry Smith   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
638729ab6d9SBarry Smith   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
639729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
640*60568592SMatthew G. Knepley   PetscStackCallHDF5(H5Oclose,(obj));
641729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(dataspace));
64236bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
64336bce6e8SMatthew G. Knepley }
64436bce6e8SMatthew G. Knepley 
645df863907SAlex Fikl /*@C
646ad0c61feSMatthew G. Knepley  PetscViewerHDF5ReadAttribute - Read a scalar attribute
647ad0c61feSMatthew G. Knepley 
648ad0c61feSMatthew G. Knepley   Input Parameters:
649ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer
650ad0c61feSMatthew G. Knepley . parent - The parent name
651ad0c61feSMatthew G. Knepley . name   - The attribute name
652ad0c61feSMatthew G. Knepley - datatype - The attribute type
653ad0c61feSMatthew G. Knepley 
654ad0c61feSMatthew G. Knepley   Output Parameter:
655ad0c61feSMatthew G. Knepley . value    - The attribute value
656ad0c61feSMatthew G. Knepley 
657ad0c61feSMatthew G. Knepley   Level: advanced
658ad0c61feSMatthew G. Knepley 
659e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute()
660ad0c61feSMatthew G. Knepley @*/
661ad0c61feSMatthew G. Knepley PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, void *value)
662ad0c61feSMatthew G. Knepley {
663*60568592SMatthew G. Knepley   hid_t          h5, obj, attribute, dtype;
664ad0c61feSMatthew G. Knepley   PetscErrorCode ierr;
665ad0c61feSMatthew G. Knepley 
666ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
6675cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
668ad0c61feSMatthew G. Knepley   PetscValidPointer(parent, 2);
669ad0c61feSMatthew G. Knepley   PetscValidPointer(name, 3);
670ad0c61feSMatthew G. Knepley   PetscValidPointer(value, 4);
671ad0c61feSMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
672ad0c61feSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
673*60568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
674*60568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
67570efba33SBarry Smith   PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
676729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
677*60568592SMatthew G. Knepley   PetscStackCallHDF5(H5Dclose,(obj));
678ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
679ad0c61feSMatthew G. Knepley }
680ad0c61feSMatthew G. Knepley 
6815cdeb108SMatthew G. Knepley static PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, const char name[], H5O_type_t otype, PetscBool *has)
6825cdeb108SMatthew G. Knepley {
6835cdeb108SMatthew G. Knepley   hid_t          h5;
6845cdeb108SMatthew G. Knepley   PetscErrorCode ierr;
6855cdeb108SMatthew G. Knepley 
6865cdeb108SMatthew G. Knepley   PetscFunctionBegin;
6875cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
6885cdeb108SMatthew G. Knepley   PetscValidPointer(name, 2);
6895cdeb108SMatthew G. Knepley   PetscValidPointer(has, 3);
6905cdeb108SMatthew G. Knepley   *has = PETSC_FALSE;
691c0cd0301SJed Brown   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
6925cdeb108SMatthew G. Knepley   if (H5Lexists(h5, name, H5P_DEFAULT)) {
6935cdeb108SMatthew G. Knepley     H5O_info_t info;
6945cdeb108SMatthew G. Knepley     hid_t      obj;
6955cdeb108SMatthew G. Knepley 
696729ab6d9SBarry Smith     PetscStackCallHDF5Return(obj,H5Oopen,(h5, name, H5P_DEFAULT));
697729ab6d9SBarry Smith     PetscStackCallHDF5(H5Oget_info,(obj, &info));
6985cdeb108SMatthew G. Knepley     if (otype == info.type) *has = PETSC_TRUE;
699729ab6d9SBarry Smith     PetscStackCallHDF5(H5Oclose,(obj));
7005cdeb108SMatthew G. Knepley   }
7015cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
7025cdeb108SMatthew G. Knepley }
7035cdeb108SMatthew G. Knepley 
704df863907SAlex Fikl /*@C
705e7bdbf8eSMatthew G. Knepley  PetscViewerHDF5HasAttribute - Check whether a scalar attribute exists
706e7bdbf8eSMatthew G. Knepley 
707e7bdbf8eSMatthew G. Knepley   Input Parameters:
708e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer
709e7bdbf8eSMatthew G. Knepley . parent - The parent name
710e7bdbf8eSMatthew G. Knepley - name   - The attribute name
711e7bdbf8eSMatthew G. Knepley 
712e7bdbf8eSMatthew G. Knepley   Output Parameter:
713e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
714e7bdbf8eSMatthew G. Knepley 
715e7bdbf8eSMatthew G. Knepley   Level: advanced
716e7bdbf8eSMatthew G. Knepley 
717e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute()
718e7bdbf8eSMatthew G. Knepley @*/
719e7bdbf8eSMatthew G. Knepley PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
720e7bdbf8eSMatthew G. Knepley {
721729ab6d9SBarry Smith   hid_t          h5, dataset;
7225cdeb108SMatthew G. Knepley   htri_t         hhas;
7235cdeb108SMatthew G. Knepley   PetscBool      exists;
724e7bdbf8eSMatthew G. Knepley   PetscErrorCode ierr;
725e7bdbf8eSMatthew G. Knepley 
726e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
7275cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
728e7bdbf8eSMatthew G. Knepley   PetscValidPointer(parent, 2);
729e7bdbf8eSMatthew G. Knepley   PetscValidPointer(name, 3);
730e7bdbf8eSMatthew G. Knepley   PetscValidPointer(has, 4);
731e7bdbf8eSMatthew G. Knepley   *has = PETSC_FALSE;
732e7bdbf8eSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
7335cdeb108SMatthew G. Knepley   ierr = PetscViewerHDF5HasObject(viewer, parent, H5O_TYPE_DATASET, &exists);CHKERRQ(ierr);
7345cdeb108SMatthew G. Knepley   if (exists) {
735e7bdbf8eSMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
736729ab6d9SBarry Smith     PetscStackCall("H5Dopen2",dataset = H5Dopen2(h5, parent, H5P_DEFAULT));
737e7bdbf8eSMatthew G. Knepley #else
738729ab6d9SBarry Smith     PetscStackCall("H5Dopen",dataset = H5Dopen(h5, parent));
739e7bdbf8eSMatthew G. Knepley #endif
740729ab6d9SBarry Smith     if (dataset < 0) PetscFunctionReturn(0);
741729ab6d9SBarry Smith     PetscStackCall("H5Aexists",hhas = H5Aexists(dataset, name));
742729ab6d9SBarry Smith     if (hhas < 0) {
743729ab6d9SBarry Smith       PetscStackCallHDF5(H5Dclose,(dataset));
744729ab6d9SBarry Smith       PetscFunctionReturn(0);
745729ab6d9SBarry Smith     }
746729ab6d9SBarry Smith     PetscStackCallHDF5(H5Dclose,(dataset));
7475cdeb108SMatthew G. Knepley     *has = hhas ? PETSC_TRUE : PETSC_FALSE;
7485cdeb108SMatthew G. Knepley   }
749e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
750e7bdbf8eSMatthew G. Knepley }
751e7bdbf8eSMatthew G. Knepley 
752a75e6a4aSMatthew G. Knepley /*
753a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
754a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
755a75e6a4aSMatthew G. Knepley */
756d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
757a75e6a4aSMatthew G. Knepley 
758a75e6a4aSMatthew G. Knepley /*@C
759a75e6a4aSMatthew G. Knepley   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
760a75e6a4aSMatthew G. Knepley 
761a75e6a4aSMatthew G. Knepley   Collective on MPI_Comm
762a75e6a4aSMatthew G. Knepley 
763a75e6a4aSMatthew G. Knepley   Input Parameter:
764a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer
765a75e6a4aSMatthew G. Knepley 
766a75e6a4aSMatthew G. Knepley   Level: intermediate
767a75e6a4aSMatthew G. Knepley 
768a75e6a4aSMatthew G. Knepley   Options Database Keys:
769a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name>
770a75e6a4aSMatthew G. Knepley 
771a75e6a4aSMatthew G. Knepley   Environmental variables:
772a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME
773a75e6a4aSMatthew G. Knepley 
774a75e6a4aSMatthew G. Knepley   Notes:
775a75e6a4aSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
776a75e6a4aSMatthew G. Knepley   an error code.  The HDF5 PetscViewer is usually used in the form
777a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
778a75e6a4aSMatthew G. Knepley 
779a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
780a75e6a4aSMatthew G. Knepley @*/
781a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
782a75e6a4aSMatthew G. Knepley {
783a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
784a75e6a4aSMatthew G. Knepley   PetscBool      flg;
785a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
786a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
787a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
788a75e6a4aSMatthew G. Knepley 
789a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
790a75e6a4aSMatthew 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);}
791a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
79212801b39SBarry Smith     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
793a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
794a75e6a4aSMatthew G. Knepley   }
79547435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
796a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
797a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
798a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
799a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
800a75e6a4aSMatthew G. Knepley     if (!flg) {
801a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname,"output.h5");
802a75e6a4aSMatthew G. Knepley       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
803a75e6a4aSMatthew G. Knepley     }
804a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
805a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
806a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
807a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
80847435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
809a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
810a75e6a4aSMatthew G. Knepley   }
811a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
812a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
813a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
814a75e6a4aSMatthew G. Knepley }
815