xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision e90c507595678ce2a232f7a34759768986da465f)
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 */
17058bd781SVaclav Hapla   char          *mataij_iname;
18058bd781SVaclav Hapla   char          *mataij_jname;
19058bd781SVaclav Hapla   char          *mataij_aname;
20058bd781SVaclav Hapla   char          *mataij_cname;
215c6c1daeSBarry Smith } PetscViewer_HDF5;
225c6c1daeSBarry Smith 
23eb5a92b4SVaclav Hapla struct _n_HDF5ReadCtx {
24eb5a92b4SVaclav Hapla   hid_t file, group, dataset, dataspace, plist;
25eb5a92b4SVaclav Hapla   PetscInt timestep;
26eb5a92b4SVaclav Hapla   PetscBool complexVal, dim2;
27eb5a92b4SVaclav Hapla };
28eb5a92b4SVaclav Hapla typedef struct _n_HDF5ReadCtx* HDF5ReadCtx;
29eb5a92b4SVaclav Hapla 
304416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v)
3182ea9b62SBarry Smith {
3282ea9b62SBarry Smith   PetscErrorCode   ierr;
3382ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;
3482ea9b62SBarry Smith 
3582ea9b62SBarry Smith   PetscFunctionBegin;
3682ea9b62SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr);
3782ea9b62SBarry Smith   ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr);
389a0502c6SHåkon Strandenes   ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr);
3982ea9b62SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
4082ea9b62SBarry Smith   PetscFunctionReturn(0);
4182ea9b62SBarry Smith }
4282ea9b62SBarry Smith 
435c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_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 = PetscFree(hdf5->filename);CHKERRQ(ierr);
50729ab6d9SBarry Smith   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
515c6c1daeSBarry Smith   PetscFunctionReturn(0);
525c6c1daeSBarry Smith }
535c6c1daeSBarry Smith 
545c6c1daeSBarry Smith PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
555c6c1daeSBarry Smith {
565c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
575c6c1daeSBarry Smith   PetscErrorCode   ierr;
585c6c1daeSBarry Smith 
595c6c1daeSBarry Smith   PetscFunctionBegin;
605c6c1daeSBarry Smith   ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr);
615c6c1daeSBarry Smith   while (hdf5->groups) {
625c6c1daeSBarry Smith     GroupList *tmp = hdf5->groups->next;
635c6c1daeSBarry Smith 
645c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups->name);CHKERRQ(ierr);
655c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups);CHKERRQ(ierr);
665c6c1daeSBarry Smith     hdf5->groups = tmp;
675c6c1daeSBarry Smith   }
685c6c1daeSBarry Smith   ierr = PetscFree(hdf5);CHKERRQ(ierr);
690b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
70d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr);
710b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
72058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr);
73058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr);
74058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetAIJNames_C",NULL);CHKERRQ(ierr);
75058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5GetAIJNames_C",NULL);CHKERRQ(ierr);
765c6c1daeSBarry Smith   PetscFunctionReturn(0);
775c6c1daeSBarry Smith }
785c6c1daeSBarry Smith 
795c6c1daeSBarry Smith PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
805c6c1daeSBarry Smith {
815c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
825c6c1daeSBarry Smith 
835c6c1daeSBarry Smith   PetscFunctionBegin;
845c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
855c6c1daeSBarry Smith   hdf5->btype = type;
865c6c1daeSBarry Smith   PetscFunctionReturn(0);
875c6c1daeSBarry Smith }
885c6c1daeSBarry Smith 
8982ea9b62SBarry Smith PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
9082ea9b62SBarry Smith {
9182ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
9282ea9b62SBarry Smith 
9382ea9b62SBarry Smith   PetscFunctionBegin;
9482ea9b62SBarry Smith   hdf5->basedimension2 = flg;
9582ea9b62SBarry Smith   PetscFunctionReturn(0);
9682ea9b62SBarry Smith }
9782ea9b62SBarry Smith 
98df863907SAlex Fikl /*@
9982ea9b62SBarry Smith      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
10082ea9b62SBarry Smith        dimension of 2.
10182ea9b62SBarry Smith 
10282ea9b62SBarry Smith     Logically Collective on PetscViewer
10382ea9b62SBarry Smith 
10482ea9b62SBarry Smith   Input Parameters:
10582ea9b62SBarry Smith +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
10682ea9b62SBarry 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
10782ea9b62SBarry Smith 
10882ea9b62SBarry Smith   Options Database:
10982ea9b62SBarry 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
11082ea9b62SBarry Smith 
11182ea9b62SBarry Smith 
11295452b02SPatrick Sanan   Notes:
11395452b02SPatrick Sanan     Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
11482ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
11582ea9b62SBarry Smith 
11682ea9b62SBarry Smith   Level: intermediate
11782ea9b62SBarry Smith 
11882ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
11982ea9b62SBarry Smith 
12082ea9b62SBarry Smith @*/
12182ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
12282ea9b62SBarry Smith {
12382ea9b62SBarry Smith   PetscErrorCode ierr;
12482ea9b62SBarry Smith 
12582ea9b62SBarry Smith   PetscFunctionBegin;
12682ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
12782ea9b62SBarry Smith   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
12882ea9b62SBarry Smith   PetscFunctionReturn(0);
12982ea9b62SBarry Smith }
13082ea9b62SBarry Smith 
131df863907SAlex Fikl /*@
13282ea9b62SBarry Smith      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
13382ea9b62SBarry Smith        dimension of 2.
13482ea9b62SBarry Smith 
13582ea9b62SBarry Smith     Logically Collective on PetscViewer
13682ea9b62SBarry Smith 
13782ea9b62SBarry Smith   Input Parameter:
13882ea9b62SBarry Smith .  viewer - the PetscViewer, must be of type HDF5
13982ea9b62SBarry Smith 
14082ea9b62SBarry Smith   Output Parameter:
14182ea9b62SBarry 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
14282ea9b62SBarry Smith 
14395452b02SPatrick Sanan   Notes:
14495452b02SPatrick Sanan     Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
14582ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
14682ea9b62SBarry Smith 
14782ea9b62SBarry Smith   Level: intermediate
14882ea9b62SBarry Smith 
14982ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
15082ea9b62SBarry Smith 
15182ea9b62SBarry Smith @*/
15282ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
15382ea9b62SBarry Smith {
15482ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
15582ea9b62SBarry Smith 
15682ea9b62SBarry Smith   PetscFunctionBegin;
15782ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
15882ea9b62SBarry Smith   *flg = hdf5->basedimension2;
15982ea9b62SBarry Smith   PetscFunctionReturn(0);
16082ea9b62SBarry Smith }
16182ea9b62SBarry Smith 
1629a0502c6SHåkon Strandenes PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
1639a0502c6SHåkon Strandenes {
1649a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1659a0502c6SHåkon Strandenes 
1669a0502c6SHåkon Strandenes   PetscFunctionBegin;
1679a0502c6SHåkon Strandenes   hdf5->spoutput = flg;
1689a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
1699a0502c6SHåkon Strandenes }
1709a0502c6SHåkon Strandenes 
171df863907SAlex Fikl /*@
1729a0502c6SHåkon Strandenes      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
1739a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
1749a0502c6SHåkon Strandenes 
1759a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
1769a0502c6SHåkon Strandenes 
1779a0502c6SHåkon Strandenes   Input Parameters:
1789a0502c6SHåkon Strandenes +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
1799a0502c6SHåkon Strandenes -  flg - if PETSC_TRUE the data will be written to disk with single precision
1809a0502c6SHåkon Strandenes 
1819a0502c6SHåkon Strandenes   Options Database:
1829a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
1839a0502c6SHåkon Strandenes 
1849a0502c6SHåkon Strandenes 
18595452b02SPatrick Sanan   Notes:
18695452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
1879a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
1889a0502c6SHåkon Strandenes 
1899a0502c6SHåkon Strandenes   Level: intermediate
1909a0502c6SHåkon Strandenes 
1919a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
1929a0502c6SHåkon Strandenes           PetscReal
1939a0502c6SHåkon Strandenes 
1949a0502c6SHåkon Strandenes @*/
1959a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
1969a0502c6SHåkon Strandenes {
1979a0502c6SHåkon Strandenes   PetscErrorCode ierr;
1989a0502c6SHåkon Strandenes 
1999a0502c6SHåkon Strandenes   PetscFunctionBegin;
2009a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2019a0502c6SHåkon Strandenes   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
2029a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2039a0502c6SHåkon Strandenes }
2049a0502c6SHåkon Strandenes 
205df863907SAlex Fikl /*@
2069a0502c6SHåkon Strandenes      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
2079a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2089a0502c6SHåkon Strandenes 
2099a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2109a0502c6SHåkon Strandenes 
2119a0502c6SHåkon Strandenes   Input Parameter:
2129a0502c6SHåkon Strandenes .  viewer - the PetscViewer, must be of type HDF5
2139a0502c6SHåkon Strandenes 
2149a0502c6SHåkon Strandenes   Output Parameter:
2159a0502c6SHåkon Strandenes .  flg - if PETSC_TRUE the data will be written to disk with single precision
2169a0502c6SHåkon Strandenes 
21795452b02SPatrick Sanan   Notes:
21895452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2199a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2209a0502c6SHåkon Strandenes 
2219a0502c6SHåkon Strandenes   Level: intermediate
2229a0502c6SHåkon Strandenes 
2239a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2249a0502c6SHåkon Strandenes           PetscReal
2259a0502c6SHåkon Strandenes 
2269a0502c6SHåkon Strandenes @*/
2279a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
2289a0502c6SHåkon Strandenes {
2299a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2309a0502c6SHåkon Strandenes 
2319a0502c6SHåkon Strandenes   PetscFunctionBegin;
2329a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2339a0502c6SHåkon Strandenes   *flg = hdf5->spoutput;
2349a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2359a0502c6SHåkon Strandenes }
2369a0502c6SHåkon Strandenes 
2375c6c1daeSBarry Smith PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
2385c6c1daeSBarry Smith {
2395c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2405c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
2415c6c1daeSBarry Smith   MPI_Info          info = MPI_INFO_NULL;
2425c6c1daeSBarry Smith #endif
2435c6c1daeSBarry Smith   hid_t             plist_id;
2445c6c1daeSBarry Smith   PetscErrorCode    ierr;
2455c6c1daeSBarry Smith 
2465c6c1daeSBarry Smith   PetscFunctionBegin;
247f683cc58SToby Isaac   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
248f683cc58SToby Isaac   if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);}
2495c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
2505c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
251729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
2525c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
253729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info));
2545c6c1daeSBarry Smith #endif
2555c6c1daeSBarry Smith   /* Create or open the file collectively */
2565c6c1daeSBarry Smith   switch (hdf5->btype) {
2575c6c1daeSBarry Smith   case FILE_MODE_READ:
258729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
2595c6c1daeSBarry Smith     break;
2605c6c1daeSBarry Smith   case FILE_MODE_APPEND:
261729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
2625c6c1daeSBarry Smith     break;
2635c6c1daeSBarry Smith   case FILE_MODE_WRITE:
264729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
2655c6c1daeSBarry Smith     break;
2665c6c1daeSBarry Smith   default:
2675c6c1daeSBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
2685c6c1daeSBarry Smith   }
2695c6c1daeSBarry Smith   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
270729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
2715c6c1daeSBarry Smith   PetscFunctionReturn(0);
2725c6c1daeSBarry Smith }
2735c6c1daeSBarry Smith 
274d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
275d1232d7fSToby Isaac {
276d1232d7fSToby Isaac   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
277d1232d7fSToby Isaac 
278d1232d7fSToby Isaac   PetscFunctionBegin;
279d1232d7fSToby Isaac   *name = vhdf5->filename;
280d1232d7fSToby Isaac   PetscFunctionReturn(0);
281d1232d7fSToby Isaac }
282d1232d7fSToby Isaac 
283058bd781SVaclav Hapla PetscErrorCode  PetscViewerHDF5SetAIJNames_HDF5(PetscViewer viewer, const char iname[], const char jname[], const char aname[], const char cname[])
284058bd781SVaclav Hapla {
285058bd781SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
286058bd781SVaclav Hapla   PetscErrorCode ierr;
287058bd781SVaclav Hapla 
288058bd781SVaclav Hapla   PetscFunctionBegin;
289058bd781SVaclav Hapla   ierr = PetscFree(hdf5->mataij_iname);CHKERRQ(ierr);
290058bd781SVaclav Hapla   ierr = PetscFree(hdf5->mataij_jname);CHKERRQ(ierr);
291058bd781SVaclav Hapla   ierr = PetscFree(hdf5->mataij_aname);CHKERRQ(ierr);
292058bd781SVaclav Hapla   ierr = PetscFree(hdf5->mataij_cname);CHKERRQ(ierr);
293058bd781SVaclav Hapla   ierr = PetscStrallocpy(iname,&hdf5->mataij_iname);CHKERRQ(ierr);
294058bd781SVaclav Hapla   ierr = PetscStrallocpy(jname,&hdf5->mataij_jname);CHKERRQ(ierr);
295058bd781SVaclav Hapla   ierr = PetscStrallocpy(aname,&hdf5->mataij_aname);CHKERRQ(ierr);
296058bd781SVaclav Hapla   ierr = PetscStrallocpy(cname,&hdf5->mataij_cname);CHKERRQ(ierr);
297058bd781SVaclav Hapla   PetscFunctionReturn(0);
298058bd781SVaclav Hapla }
299058bd781SVaclav Hapla 
300058bd781SVaclav Hapla /*@C
301058bd781SVaclav Hapla   PetscViewerHDF5SetAIJNames - Set the names of the datasets representing the three AIJ (CRS) arrays and the name of the attribute storing the number of columns within the HDF5 file.
302058bd781SVaclav Hapla 
303058bd781SVaclav Hapla   Collective on PetscViewer
304058bd781SVaclav Hapla 
305058bd781SVaclav Hapla   Input Parameters:
306058bd781SVaclav Hapla +  viewer - the PetscViewer; either ASCII or binary
307058bd781SVaclav Hapla .  iname - name of dataset i representing row pointers; that is i[0] = 0, i[row] = i[row-1] + number of elements in that row of the matrix
308058bd781SVaclav Hapla .  jname - name of dataset j representing column indices
309058bd781SVaclav Hapla .  aname - name of dataset a representing matrix values
310058bd781SVaclav Hapla -  cname - name of attribute stoting column count
311058bd781SVaclav Hapla 
312058bd781SVaclav Hapla   Level: advanced
313058bd781SVaclav Hapla 
314058bd781SVaclav Hapla   Notes:
315058bd781SVaclav Hapla   Current defaults are (iname, jname, aname, cname) = ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be readily loaded.
316058bd781SVaclav Hapla 
317058bd781SVaclav Hapla .seealso: MatLoad(), PetscViewerCreate(), PetscViewerSetType(), PETSCVIEWERHDF5, PetscViewerHDF5GetAIJNames()
318058bd781SVaclav Hapla @*/
319058bd781SVaclav Hapla PetscErrorCode  PetscViewerHDF5SetAIJNames(PetscViewer viewer, const char iname[], const char jname[], const char aname[], const char cname[])
320058bd781SVaclav Hapla {
321058bd781SVaclav Hapla   PetscErrorCode ierr;
322058bd781SVaclav Hapla 
323058bd781SVaclav Hapla   PetscFunctionBegin;
324058bd781SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
325058bd781SVaclav Hapla   PetscValidCharPointer(iname,2);
326058bd781SVaclav Hapla   PetscValidCharPointer(jname,3);
327058bd781SVaclav Hapla   PetscValidCharPointer(aname,4);
328058bd781SVaclav Hapla   PetscValidCharPointer(cname,5);
329058bd781SVaclav Hapla   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetAIJNames_C",(PetscViewer,const char[],const char[],const char[],const char[]),(viewer,iname,jname,aname,cname));CHKERRQ(ierr);
330058bd781SVaclav Hapla   PetscFunctionReturn(0);
331058bd781SVaclav Hapla }
332058bd781SVaclav Hapla 
333058bd781SVaclav Hapla PetscErrorCode  PetscViewerHDF5GetAIJNames_HDF5(PetscViewer viewer, const char *iname[], const char *jname[], const char *aname[], const char *cname[])
334058bd781SVaclav Hapla {
335058bd781SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
336058bd781SVaclav Hapla 
337058bd781SVaclav Hapla   PetscFunctionBegin;
338058bd781SVaclav Hapla   *iname = hdf5->mataij_iname;
339058bd781SVaclav Hapla   *jname = hdf5->mataij_jname;
340058bd781SVaclav Hapla   *aname = hdf5->mataij_aname;
341058bd781SVaclav Hapla   *cname = hdf5->mataij_cname;
342058bd781SVaclav Hapla   PetscFunctionReturn(0);
343058bd781SVaclav Hapla }
344058bd781SVaclav Hapla 
345058bd781SVaclav Hapla /*@C
346058bd781SVaclav Hapla   PetscViewerHDF5GetAIJNames - Get the names of the datasets representing the three AIJ (CRS) arrays and the name of the attribute storing the number of columns within the HDF5 file.
347058bd781SVaclav Hapla 
348058bd781SVaclav Hapla   Collective on PetscViewer
349058bd781SVaclav Hapla 
350058bd781SVaclav Hapla   Input Parameters:
351058bd781SVaclav Hapla .  viewer - the PetscViewer; either ASCII or binary
352058bd781SVaclav Hapla 
353058bd781SVaclav Hapla   Output Parameters:
354058bd781SVaclav Hapla +  iname - name of dataset i representing row pointers; that is i[0] = 0, i[row] = i[row-1] + number of elements in that row of the matrix
355058bd781SVaclav Hapla .  jname - name of dataset j representing column indices
356058bd781SVaclav Hapla .  aname - name of dataset a representing matrix values
357058bd781SVaclav Hapla -  cname - name of attribute stoting column count
358058bd781SVaclav Hapla 
359058bd781SVaclav Hapla   Level: advanced
360058bd781SVaclav Hapla 
361058bd781SVaclav Hapla   Notes:
362058bd781SVaclav Hapla   Current defaults are (iname, jname, aname, cname) = ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be readily loaded.
363058bd781SVaclav Hapla 
364058bd781SVaclav Hapla .seealso: MatLoad(), PetscViewerCreate(), PetscViewerSetType(), PETSCVIEWERHDF5, PetscViewerHDF5SetAIJNames()
365058bd781SVaclav Hapla @*/
366058bd781SVaclav Hapla PetscErrorCode  PetscViewerHDF5GetAIJNames(PetscViewer viewer, const char *iname[], const char *jname[], const char *aname[], const char *cname[])
367058bd781SVaclav Hapla {
368058bd781SVaclav Hapla   PetscErrorCode ierr;
369058bd781SVaclav Hapla 
370058bd781SVaclav Hapla   PetscFunctionBegin;
371058bd781SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
372058bd781SVaclav Hapla   PetscValidPointer(iname,2);
373058bd781SVaclav Hapla   PetscValidPointer(jname,3);
374058bd781SVaclav Hapla   PetscValidPointer(aname,4);
375058bd781SVaclav Hapla   PetscValidPointer(cname,5);
376058bd781SVaclav Hapla   ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetAIJNames_C",(PetscViewer,const char*[],const char*[],const char*[],const char*[]),(viewer,iname,jname,aname,cname));CHKERRQ(ierr);
377058bd781SVaclav Hapla   PetscFunctionReturn(0);
378058bd781SVaclav Hapla }
379058bd781SVaclav Hapla 
3808556b5ebSBarry Smith /*MC
3818556b5ebSBarry Smith    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
3828556b5ebSBarry Smith 
3838556b5ebSBarry Smith 
3848556b5ebSBarry Smith .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
3858556b5ebSBarry Smith            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
3868556b5ebSBarry Smith            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
3878556b5ebSBarry Smith            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
3888556b5ebSBarry Smith 
3891b266c99SBarry Smith   Level: beginner
3908556b5ebSBarry Smith M*/
391d1232d7fSToby Isaac 
3928cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
3935c6c1daeSBarry Smith {
3945c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
3955c6c1daeSBarry Smith   PetscErrorCode   ierr;
3965c6c1daeSBarry Smith 
3975c6c1daeSBarry Smith   PetscFunctionBegin;
398b00a9115SJed Brown   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
3995c6c1daeSBarry Smith 
4005c6c1daeSBarry Smith   v->data                = (void*) hdf5;
4015c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
40282ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
4035c6c1daeSBarry Smith   v->ops->flush          = 0;
4045c6c1daeSBarry Smith   hdf5->btype            = (PetscFileMode) -1;
4055c6c1daeSBarry Smith   hdf5->filename         = 0;
4065c6c1daeSBarry Smith   hdf5->timestep         = -1;
4070298fd71SBarry Smith   hdf5->groups           = NULL;
4085c6c1daeSBarry Smith 
409058bd781SVaclav Hapla   /* ir and jc are deliberately swapped as MATLAB uses column-major format */
410058bd781SVaclav Hapla   ierr = PetscStrallocpy("jc",  &hdf5->mataij_iname);CHKERRQ(ierr);
411058bd781SVaclav Hapla   ierr = PetscStrallocpy("ir",  &hdf5->mataij_jname);CHKERRQ(ierr);
412058bd781SVaclav Hapla   ierr = PetscStrallocpy("data",&hdf5->mataij_aname);CHKERRQ(ierr);
413058bd781SVaclav Hapla   ierr = PetscStrallocpy("MATLAB_sparse", &hdf5->mataij_cname);CHKERRQ(ierr);
414058bd781SVaclav Hapla 
4150b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
416d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr);
4170b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
41882ea9b62SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
4199a0502c6SHåkon Strandenes   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr);
420058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetAIJNames_C",PetscViewerHDF5SetAIJNames_HDF5);CHKERRQ(ierr);
421058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetAIJNames_C",PetscViewerHDF5GetAIJNames_HDF5);CHKERRQ(ierr);
4225c6c1daeSBarry Smith   PetscFunctionReturn(0);
4235c6c1daeSBarry Smith }
4245c6c1daeSBarry Smith 
4255c6c1daeSBarry Smith /*@C
4265c6c1daeSBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
4275c6c1daeSBarry Smith 
4285c6c1daeSBarry Smith    Collective on MPI_Comm
4295c6c1daeSBarry Smith 
4305c6c1daeSBarry Smith    Input Parameters:
4315c6c1daeSBarry Smith +  comm - MPI communicator
4325c6c1daeSBarry Smith .  name - name of file
4335c6c1daeSBarry Smith -  type - type of file
4345c6c1daeSBarry Smith $    FILE_MODE_WRITE - create new file for binary output
4355c6c1daeSBarry Smith $    FILE_MODE_READ - open existing file for binary input
4365c6c1daeSBarry Smith $    FILE_MODE_APPEND - open existing file for binary output
4375c6c1daeSBarry Smith 
4385c6c1daeSBarry Smith    Output Parameter:
4395c6c1daeSBarry Smith .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
4405c6c1daeSBarry Smith 
44182ea9b62SBarry Smith   Options Database:
44282ea9b62SBarry 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
4439a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
44482ea9b62SBarry Smith 
4455c6c1daeSBarry Smith    Level: beginner
4465c6c1daeSBarry Smith 
4475c6c1daeSBarry Smith    Note:
4485c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
4495c6c1daeSBarry Smith 
4505c6c1daeSBarry Smith    Concepts: HDF5 files
4515c6c1daeSBarry Smith    Concepts: PetscViewerHDF5^creating
4525c6c1daeSBarry Smith 
4536a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
4549a0502c6SHåkon Strandenes           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
455a56f64adSBarry Smith           MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
4565c6c1daeSBarry Smith @*/
4575c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
4585c6c1daeSBarry Smith {
4595c6c1daeSBarry Smith   PetscErrorCode ierr;
4605c6c1daeSBarry Smith 
4615c6c1daeSBarry Smith   PetscFunctionBegin;
4625c6c1daeSBarry Smith   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
4635c6c1daeSBarry Smith   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
4645c6c1daeSBarry Smith   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
4655c6c1daeSBarry Smith   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
4665c6c1daeSBarry Smith   PetscFunctionReturn(0);
4675c6c1daeSBarry Smith }
4685c6c1daeSBarry Smith 
4695c6c1daeSBarry Smith /*@C
4705c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
4715c6c1daeSBarry Smith 
4725c6c1daeSBarry Smith   Not collective
4735c6c1daeSBarry Smith 
4745c6c1daeSBarry Smith   Input Parameter:
4755c6c1daeSBarry Smith . viewer - the PetscViewer
4765c6c1daeSBarry Smith 
4775c6c1daeSBarry Smith   Output Parameter:
4785c6c1daeSBarry Smith . file_id - The file id
4795c6c1daeSBarry Smith 
4805c6c1daeSBarry Smith   Level: intermediate
4815c6c1daeSBarry Smith 
4825c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open()
4835c6c1daeSBarry Smith @*/
4845c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
4855c6c1daeSBarry Smith {
4865c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4875c6c1daeSBarry Smith 
4885c6c1daeSBarry Smith   PetscFunctionBegin;
4895c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
4905c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
4915c6c1daeSBarry Smith   PetscFunctionReturn(0);
4925c6c1daeSBarry Smith }
4935c6c1daeSBarry Smith 
4945c6c1daeSBarry Smith /*@C
4955c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
4965c6c1daeSBarry Smith 
4975c6c1daeSBarry Smith   Not collective
4985c6c1daeSBarry Smith 
4995c6c1daeSBarry Smith   Input Parameters:
5005c6c1daeSBarry Smith + viewer - the PetscViewer
5015c6c1daeSBarry Smith - name - The group name
5025c6c1daeSBarry Smith 
5035c6c1daeSBarry Smith   Level: intermediate
5045c6c1daeSBarry Smith 
505874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
5065c6c1daeSBarry Smith @*/
5075c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
5085c6c1daeSBarry Smith {
5095c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5105c6c1daeSBarry Smith   GroupList        *groupNode;
5115c6c1daeSBarry Smith   PetscErrorCode   ierr;
5125c6c1daeSBarry Smith 
5135c6c1daeSBarry Smith   PetscFunctionBegin;
5145c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5155c6c1daeSBarry Smith   PetscValidCharPointer(name,2);
51695dccacaSBarry Smith   ierr = PetscNew(&groupNode);CHKERRQ(ierr);
5175c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
518a297a907SKarl Rupp 
5195c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
5205c6c1daeSBarry Smith   hdf5->groups    = groupNode;
5215c6c1daeSBarry Smith   PetscFunctionReturn(0);
5225c6c1daeSBarry Smith }
5235c6c1daeSBarry Smith 
5243ef9c667SSatish Balay /*@
5255c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
5265c6c1daeSBarry Smith 
5275c6c1daeSBarry Smith   Not collective
5285c6c1daeSBarry Smith 
5295c6c1daeSBarry Smith   Input Parameter:
5305c6c1daeSBarry Smith . viewer - the PetscViewer
5315c6c1daeSBarry Smith 
5325c6c1daeSBarry Smith   Level: intermediate
5335c6c1daeSBarry Smith 
534874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
5355c6c1daeSBarry Smith @*/
5365c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
5375c6c1daeSBarry Smith {
5385c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5395c6c1daeSBarry Smith   GroupList        *groupNode;
5405c6c1daeSBarry Smith   PetscErrorCode   ierr;
5415c6c1daeSBarry Smith 
5425c6c1daeSBarry Smith   PetscFunctionBegin;
5435c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
54482f516ccSBarry Smith   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
5455c6c1daeSBarry Smith   groupNode    = hdf5->groups;
5465c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
5475c6c1daeSBarry Smith   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
5485c6c1daeSBarry Smith   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
5495c6c1daeSBarry Smith   PetscFunctionReturn(0);
5505c6c1daeSBarry Smith }
5515c6c1daeSBarry Smith 
5525c6c1daeSBarry Smith /*@C
553874270d9SVaclav Hapla   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup().
554874270d9SVaclav Hapla   If none has been assigned, returns NULL.
5555c6c1daeSBarry Smith 
5565c6c1daeSBarry Smith   Not collective
5575c6c1daeSBarry Smith 
5585c6c1daeSBarry Smith   Input Parameter:
5595c6c1daeSBarry Smith . viewer - the PetscViewer
5605c6c1daeSBarry Smith 
5615c6c1daeSBarry Smith   Output Parameter:
5625c6c1daeSBarry Smith . name - The group name
5635c6c1daeSBarry Smith 
5645c6c1daeSBarry Smith   Level: intermediate
5655c6c1daeSBarry Smith 
566874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
5675c6c1daeSBarry Smith @*/
5685c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
5695c6c1daeSBarry Smith {
5705c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
5715c6c1daeSBarry Smith 
5725c6c1daeSBarry Smith   PetscFunctionBegin;
5735c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
574c959eef4SJed Brown   PetscValidPointer(name,2);
575a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
5760298fd71SBarry Smith   else *name = NULL;
5775c6c1daeSBarry Smith   PetscFunctionReturn(0);
5785c6c1daeSBarry Smith }
5795c6c1daeSBarry Smith 
5808c557b5aSVaclav Hapla /*@
581874270d9SVaclav Hapla   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(),
582874270d9SVaclav Hapla   and return this group's ID and file ID.
583874270d9SVaclav Hapla   If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID.
584874270d9SVaclav Hapla 
585874270d9SVaclav Hapla   Not collective
586874270d9SVaclav Hapla 
587874270d9SVaclav Hapla   Input Parameter:
588874270d9SVaclav Hapla . viewer - the PetscViewer
589874270d9SVaclav Hapla 
590874270d9SVaclav Hapla   Output Parameter:
591874270d9SVaclav Hapla + fileId - The HDF5 file ID
592874270d9SVaclav Hapla - groupId - The HDF5 group ID
593874270d9SVaclav Hapla 
594874270d9SVaclav Hapla   Level: intermediate
595874270d9SVaclav Hapla 
596874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
597874270d9SVaclav Hapla @*/
59854dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
59954dbf706SVaclav Hapla {
60054dbf706SVaclav Hapla   hid_t          file_id, group;
60154dbf706SVaclav Hapla   htri_t         found;
60254dbf706SVaclav Hapla   const char     *groupName = NULL;
60354dbf706SVaclav Hapla   PetscErrorCode ierr;
60454dbf706SVaclav Hapla 
60554dbf706SVaclav Hapla   PetscFunctionBegin;
60654dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr);
60754dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr);
60854dbf706SVaclav Hapla   /* Open group */
60954dbf706SVaclav Hapla   if (groupName) {
61054dbf706SVaclav Hapla     PetscBool root;
61154dbf706SVaclav Hapla 
61254dbf706SVaclav Hapla     ierr = PetscStrcmp(groupName, "/", &root);CHKERRQ(ierr);
61354dbf706SVaclav Hapla     PetscStackCall("H5Lexists",found = H5Lexists(file_id, groupName, H5P_DEFAULT));
61454dbf706SVaclav Hapla     if (!root && (found <= 0)) {
61554dbf706SVaclav Hapla #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
61654dbf706SVaclav Hapla       PetscStackCallHDF5Return(group,H5Gcreate2,(file_id, groupName, 0, H5P_DEFAULT, H5P_DEFAULT));
61754dbf706SVaclav Hapla #else /* deprecated HDF5 1.6 API */
61854dbf706SVaclav Hapla       PetscStackCallHDF5Return(group,H5Gcreate,(file_id, groupName, 0));
61954dbf706SVaclav Hapla #endif
62054dbf706SVaclav Hapla       PetscStackCallHDF5(H5Gclose,(group));
62154dbf706SVaclav Hapla     }
62254dbf706SVaclav Hapla #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
62354dbf706SVaclav Hapla     PetscStackCallHDF5Return(group,H5Gopen2,(file_id, groupName, H5P_DEFAULT));
62454dbf706SVaclav Hapla #else
62554dbf706SVaclav Hapla     PetscStackCallHDF5Return(group,H5Gopen,(file_id, groupName));
62654dbf706SVaclav Hapla #endif
62754dbf706SVaclav Hapla   } else group = file_id;
62854dbf706SVaclav Hapla 
62954dbf706SVaclav Hapla   *fileId  = file_id;
63054dbf706SVaclav Hapla   *groupId = group;
63154dbf706SVaclav Hapla   PetscFunctionReturn(0);
63254dbf706SVaclav Hapla }
63354dbf706SVaclav Hapla 
6343ef9c667SSatish Balay /*@
6355c6c1daeSBarry Smith   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
6365c6c1daeSBarry Smith 
6375c6c1daeSBarry Smith   Not collective
6385c6c1daeSBarry Smith 
6395c6c1daeSBarry Smith   Input Parameter:
6405c6c1daeSBarry Smith . viewer - the PetscViewer
6415c6c1daeSBarry Smith 
6425c6c1daeSBarry Smith   Level: intermediate
6435c6c1daeSBarry Smith 
6445c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
6455c6c1daeSBarry Smith @*/
6465c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
6475c6c1daeSBarry Smith {
6485c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
6495c6c1daeSBarry Smith 
6505c6c1daeSBarry Smith   PetscFunctionBegin;
6515c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
6525c6c1daeSBarry Smith   ++hdf5->timestep;
6535c6c1daeSBarry Smith   PetscFunctionReturn(0);
6545c6c1daeSBarry Smith }
6555c6c1daeSBarry Smith 
6563ef9c667SSatish Balay /*@
6575c6c1daeSBarry Smith   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
6585c6c1daeSBarry Smith   of -1 disables blocking with timesteps.
6595c6c1daeSBarry Smith 
6605c6c1daeSBarry Smith   Not collective
6615c6c1daeSBarry Smith 
6625c6c1daeSBarry Smith   Input Parameters:
6635c6c1daeSBarry Smith + viewer - the PetscViewer
6645c6c1daeSBarry Smith - timestep - The timestep number
6655c6c1daeSBarry Smith 
6665c6c1daeSBarry Smith   Level: intermediate
6675c6c1daeSBarry Smith 
6685c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
6695c6c1daeSBarry Smith @*/
6705c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
6715c6c1daeSBarry Smith {
6725c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
6735c6c1daeSBarry Smith 
6745c6c1daeSBarry Smith   PetscFunctionBegin;
6755c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
6765c6c1daeSBarry Smith   hdf5->timestep = timestep;
6775c6c1daeSBarry Smith   PetscFunctionReturn(0);
6785c6c1daeSBarry Smith }
6795c6c1daeSBarry Smith 
6803ef9c667SSatish Balay /*@
6815c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
6825c6c1daeSBarry Smith 
6835c6c1daeSBarry Smith   Not collective
6845c6c1daeSBarry Smith 
6855c6c1daeSBarry Smith   Input Parameter:
6865c6c1daeSBarry Smith . viewer - the PetscViewer
6875c6c1daeSBarry Smith 
6885c6c1daeSBarry Smith   Output Parameter:
6895c6c1daeSBarry Smith . timestep - The timestep number
6905c6c1daeSBarry Smith 
6915c6c1daeSBarry Smith   Level: intermediate
6925c6c1daeSBarry Smith 
6935c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
6945c6c1daeSBarry Smith @*/
6955c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
6965c6c1daeSBarry Smith {
6975c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
6985c6c1daeSBarry Smith 
6995c6c1daeSBarry Smith   PetscFunctionBegin;
7005c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7015c6c1daeSBarry Smith   PetscValidPointer(timestep,2);
7025c6c1daeSBarry Smith   *timestep = hdf5->timestep;
7035c6c1daeSBarry Smith   PetscFunctionReturn(0);
7045c6c1daeSBarry Smith }
7055c6c1daeSBarry Smith 
70636bce6e8SMatthew G. Knepley /*@C
70736bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
70836bce6e8SMatthew G. Knepley 
70936bce6e8SMatthew G. Knepley   Not collective
71036bce6e8SMatthew G. Knepley 
71136bce6e8SMatthew G. Knepley   Input Parameter:
71236bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
71336bce6e8SMatthew G. Knepley 
71436bce6e8SMatthew G. Knepley   Output Parameter:
71536bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
71636bce6e8SMatthew G. Knepley 
71736bce6e8SMatthew G. Knepley   Level: advanced
71836bce6e8SMatthew G. Knepley 
71936bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
72036bce6e8SMatthew G. Knepley @*/
72136bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
72236bce6e8SMatthew G. Knepley {
72336bce6e8SMatthew G. Knepley   PetscFunctionBegin;
72436bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
72536bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
72636bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_LLONG;
72736bce6e8SMatthew G. Knepley #else
72836bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_INT;
72936bce6e8SMatthew G. Knepley #endif
73036bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
73136bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
73236bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
73336bce6e8SMatthew G. Knepley   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_DOUBLE;
73436bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_DOUBLE;
73536bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
73636bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
73736bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
7387e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
73936bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
74036bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
74136bce6e8SMatthew G. Knepley }
74236bce6e8SMatthew G. Knepley 
74336bce6e8SMatthew G. Knepley /*@C
74436bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
74536bce6e8SMatthew G. Knepley 
74636bce6e8SMatthew G. Knepley   Not collective
74736bce6e8SMatthew G. Knepley 
74836bce6e8SMatthew G. Knepley   Input Parameter:
74936bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
75036bce6e8SMatthew G. Knepley 
75136bce6e8SMatthew G. Knepley   Output Parameter:
75236bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
75336bce6e8SMatthew G. Knepley 
75436bce6e8SMatthew G. Knepley   Level: advanced
75536bce6e8SMatthew G. Knepley 
75636bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
75736bce6e8SMatthew G. Knepley @*/
75836bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
75936bce6e8SMatthew G. Knepley {
76036bce6e8SMatthew G. Knepley   PetscFunctionBegin;
76136bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
76236bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
76336bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
76436bce6e8SMatthew G. Knepley #else
76536bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
76636bce6e8SMatthew G. Knepley #endif
76736bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
76836bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
76936bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
77036bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
77136bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
77236bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
7737e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
77436bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
77536bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
77636bce6e8SMatthew G. Knepley }
77736bce6e8SMatthew G. Knepley 
778df863907SAlex Fikl /*@C
77936bce6e8SMatthew G. Knepley  PetscViewerHDF5WriteAttribute - Write a scalar attribute
78036bce6e8SMatthew G. Knepley 
78136bce6e8SMatthew G. Knepley   Input Parameters:
78236bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer
78336bce6e8SMatthew G. Knepley . parent - The parent name
78436bce6e8SMatthew G. Knepley . name   - The attribute name
78536bce6e8SMatthew G. Knepley . datatype - The attribute type
78636bce6e8SMatthew G. Knepley - value    - The attribute value
78736bce6e8SMatthew G. Knepley 
78836bce6e8SMatthew G. Knepley   Level: advanced
78936bce6e8SMatthew G. Knepley 
790e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute()
79136bce6e8SMatthew G. Knepley @*/
79236bce6e8SMatthew G. Knepley PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
79336bce6e8SMatthew G. Knepley {
79460568592SMatthew G. Knepley   hid_t          h5, dataspace, obj, attribute, dtype;
79536bce6e8SMatthew G. Knepley   PetscErrorCode ierr;
79636bce6e8SMatthew G. Knepley 
79736bce6e8SMatthew G. Knepley   PetscFunctionBegin;
7985cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
79936bce6e8SMatthew G. Knepley   PetscValidPointer(parent, 2);
80036bce6e8SMatthew G. Knepley   PetscValidPointer(name, 3);
80136bce6e8SMatthew G. Knepley   PetscValidPointer(value, 4);
80236bce6e8SMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
8037e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
8047e97c476SMatthew G. Knepley     size_t len;
8057e97c476SMatthew G. Knepley     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
806729ab6d9SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
8077e97c476SMatthew G. Knepley   }
80836bce6e8SMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
809729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
81060568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
81136bce6e8SMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
81260568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
81336bce6e8SMatthew G. Knepley #else
81460568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Acreate,(obj, name, dtype, dataspace, H5P_DEFAULT));
81536bce6e8SMatthew G. Knepley #endif
816729ab6d9SBarry Smith   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
817729ab6d9SBarry Smith   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
818729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
81960568592SMatthew G. Knepley   PetscStackCallHDF5(H5Oclose,(obj));
820729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(dataspace));
82136bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
82236bce6e8SMatthew G. Knepley }
82336bce6e8SMatthew G. Knepley 
824df863907SAlex Fikl /*@C
825ad0c61feSMatthew G. Knepley  PetscViewerHDF5ReadAttribute - Read a scalar attribute
826ad0c61feSMatthew G. Knepley 
827ad0c61feSMatthew G. Knepley   Input Parameters:
828ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer
829ad0c61feSMatthew G. Knepley . parent - The parent name
830ad0c61feSMatthew G. Knepley . name   - The attribute name
831ad0c61feSMatthew G. Knepley - datatype - The attribute type
832ad0c61feSMatthew G. Knepley 
833ad0c61feSMatthew G. Knepley   Output Parameter:
834ad0c61feSMatthew G. Knepley . value    - The attribute value
835ad0c61feSMatthew G. Knepley 
836ad0c61feSMatthew G. Knepley   Level: advanced
837ad0c61feSMatthew G. Knepley 
838e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute()
839ad0c61feSMatthew G. Knepley @*/
840ad0c61feSMatthew G. Knepley PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, void *value)
841ad0c61feSMatthew G. Knepley {
842f0b58479SMatthew G. Knepley   hid_t          h5, obj, attribute, atype, dtype;
843ad0c61feSMatthew G. Knepley   PetscErrorCode ierr;
844ad0c61feSMatthew G. Knepley 
845ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
8465cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
847ad0c61feSMatthew G. Knepley   PetscValidPointer(parent, 2);
848ad0c61feSMatthew G. Knepley   PetscValidPointer(name, 3);
849ad0c61feSMatthew G. Knepley   PetscValidPointer(value, 4);
850ad0c61feSMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
851ad0c61feSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
85260568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
85360568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
854f0b58479SMatthew G. Knepley   if (datatype == PETSC_STRING) {
855f0b58479SMatthew G. Knepley     size_t len;
856*e90c5075SMarek Pecha     PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
857f0b58479SMatthew G. Knepley     PetscStackCallHDF5Return(len,H5Tget_size,(atype));
858f0b58479SMatthew G. Knepley     PetscStackCallHDF5(H5Tclose,(atype));
859f0b58479SMatthew G. Knepley     ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr);
860f0b58479SMatthew G. Knepley   }
86170efba33SBarry Smith   PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
862729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
863*e90c5075SMarek Pecha   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
864*e90c5075SMarek Pecha   PetscStackCallHDF5(H5Oclose,(obj));
865ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
866ad0c61feSMatthew G. Knepley }
867ad0c61feSMatthew G. Knepley 
8685cdeb108SMatthew G. Knepley static PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, const char name[], H5O_type_t otype, PetscBool *has)
8695cdeb108SMatthew G. Knepley {
8705cdeb108SMatthew G. Knepley   hid_t          h5;
8715cdeb108SMatthew G. Knepley   PetscErrorCode ierr;
8725cdeb108SMatthew G. Knepley 
8735cdeb108SMatthew G. Knepley   PetscFunctionBegin;
8745cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
8755cdeb108SMatthew G. Knepley   PetscValidPointer(name, 2);
8765cdeb108SMatthew G. Knepley   PetscValidPointer(has, 3);
8775cdeb108SMatthew G. Knepley   *has = PETSC_FALSE;
878c0cd0301SJed Brown   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
8795cdeb108SMatthew G. Knepley   if (H5Lexists(h5, name, H5P_DEFAULT)) {
8805cdeb108SMatthew G. Knepley     H5O_info_t info;
8815cdeb108SMatthew G. Knepley     hid_t      obj;
8825cdeb108SMatthew G. Knepley 
883729ab6d9SBarry Smith     PetscStackCallHDF5Return(obj,H5Oopen,(h5, name, H5P_DEFAULT));
884729ab6d9SBarry Smith     PetscStackCallHDF5(H5Oget_info,(obj, &info));
8855cdeb108SMatthew G. Knepley     if (otype == info.type) *has = PETSC_TRUE;
886729ab6d9SBarry Smith     PetscStackCallHDF5(H5Oclose,(obj));
8875cdeb108SMatthew G. Knepley   }
8885cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
8895cdeb108SMatthew G. Knepley }
8905cdeb108SMatthew G. Knepley 
891df863907SAlex Fikl /*@C
892e7bdbf8eSMatthew G. Knepley  PetscViewerHDF5HasAttribute - Check whether a scalar attribute exists
893e7bdbf8eSMatthew G. Knepley 
894e7bdbf8eSMatthew G. Knepley   Input Parameters:
895e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer
896e7bdbf8eSMatthew G. Knepley . parent - The parent name
897e7bdbf8eSMatthew G. Knepley - name   - The attribute name
898e7bdbf8eSMatthew G. Knepley 
899e7bdbf8eSMatthew G. Knepley   Output Parameter:
900e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
901e7bdbf8eSMatthew G. Knepley 
902e7bdbf8eSMatthew G. Knepley   Level: advanced
903e7bdbf8eSMatthew G. Knepley 
904e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute()
905e7bdbf8eSMatthew G. Knepley @*/
906e7bdbf8eSMatthew G. Knepley PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
907e7bdbf8eSMatthew G. Knepley {
908729ab6d9SBarry Smith   hid_t          h5, dataset;
9095cdeb108SMatthew G. Knepley   htri_t         hhas;
9105cdeb108SMatthew G. Knepley   PetscBool      exists;
911e7bdbf8eSMatthew G. Knepley   PetscErrorCode ierr;
912e7bdbf8eSMatthew G. Knepley 
913e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
9145cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
915e7bdbf8eSMatthew G. Knepley   PetscValidPointer(parent, 2);
916e7bdbf8eSMatthew G. Knepley   PetscValidPointer(name, 3);
917e7bdbf8eSMatthew G. Knepley   PetscValidPointer(has, 4);
918e7bdbf8eSMatthew G. Knepley   *has = PETSC_FALSE;
919e7bdbf8eSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
9205cdeb108SMatthew G. Knepley   ierr = PetscViewerHDF5HasObject(viewer, parent, H5O_TYPE_DATASET, &exists);CHKERRQ(ierr);
9215cdeb108SMatthew G. Knepley   if (exists) {
922e7bdbf8eSMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
923729ab6d9SBarry Smith     PetscStackCall("H5Dopen2",dataset = H5Dopen2(h5, parent, H5P_DEFAULT));
924e7bdbf8eSMatthew G. Knepley #else
925729ab6d9SBarry Smith     PetscStackCall("H5Dopen",dataset = H5Dopen(h5, parent));
926e7bdbf8eSMatthew G. Knepley #endif
927729ab6d9SBarry Smith     if (dataset < 0) PetscFunctionReturn(0);
928729ab6d9SBarry Smith     PetscStackCall("H5Aexists",hhas = H5Aexists(dataset, name));
929729ab6d9SBarry Smith     if (hhas < 0) {
930729ab6d9SBarry Smith       PetscStackCallHDF5(H5Dclose,(dataset));
931729ab6d9SBarry Smith       PetscFunctionReturn(0);
932729ab6d9SBarry Smith     }
933729ab6d9SBarry Smith     PetscStackCallHDF5(H5Dclose,(dataset));
9345cdeb108SMatthew G. Knepley     *has = hhas ? PETSC_TRUE : PETSC_FALSE;
9355cdeb108SMatthew G. Knepley   }
936e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
937e7bdbf8eSMatthew G. Knepley }
938e7bdbf8eSMatthew G. Knepley 
939eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadInitialize_Private(PetscViewer viewer, const char name[], HDF5ReadCtx *ctx)
940a5e1feadSVaclav Hapla {
941fbd37863SVaclav Hapla   HDF5ReadCtx    h=NULL;
942fbd37863SVaclav Hapla   const char    *groupname=NULL;
94369a06e7bSVaclav Hapla   char           vecgroup[PETSC_MAX_PATH_LEN];
944a5e1feadSVaclav Hapla   PetscErrorCode ierr;
945a5e1feadSVaclav Hapla 
946a5e1feadSVaclav Hapla   PetscFunctionBegin;
94769a06e7bSVaclav Hapla   ierr = PetscNew(&h);CHKERRQ(ierr);
94869a06e7bSVaclav Hapla   ierr = PetscViewerHDF5OpenGroup(viewer, &h->file, &h->group);CHKERRQ(ierr);
949a5e1feadSVaclav Hapla #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
95069a06e7bSVaclav Hapla   PetscStackCallHDF5Return(h->dataset,H5Dopen2,(h->group, name, H5P_DEFAULT));
951a5e1feadSVaclav Hapla #else
95269a06e7bSVaclav Hapla   PetscStackCallHDF5Return(h->dataset,H5Dopen,(h->group, name));
953a5e1feadSVaclav Hapla #endif
95469a06e7bSVaclav Hapla   PetscStackCallHDF5Return(h->dataspace,H5Dget_space,(h->dataset));
9559568d583SVaclav Hapla   ierr = PetscViewerHDF5GetTimestep(viewer, &h->timestep);CHKERRQ(ierr);
95669a06e7bSVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer,&groupname);CHKERRQ(ierr);
957b44ade6dSVaclav Hapla   ierr = PetscSNPrintf(vecgroup,PETSC_MAX_PATH_LEN,"%s/%s",groupname ? groupname : "",name);CHKERRQ(ierr);
9589568d583SVaclav Hapla   ierr = PetscViewerHDF5HasAttribute(viewer,vecgroup,"complex",&h->complexVal);CHKERRQ(ierr);
959b8ef406cSVaclav Hapla   /* Create property list for collective dataset read */
960b8ef406cSVaclav Hapla   PetscStackCallHDF5Return(h->plist,H5Pcreate,(H5P_DATASET_XFER));
961b8ef406cSVaclav Hapla #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
962b8ef406cSVaclav Hapla   PetscStackCallHDF5(H5Pset_dxpl_mpio,(h->plist, H5FD_MPIO_COLLECTIVE));
963b8ef406cSVaclav Hapla #endif
964b8ef406cSVaclav Hapla   /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */
96569a06e7bSVaclav Hapla   *ctx = h;
96669a06e7bSVaclav Hapla   PetscFunctionReturn(0);
96769a06e7bSVaclav Hapla }
96869a06e7bSVaclav Hapla 
969eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadFinalize_Private(PetscViewer viewer, HDF5ReadCtx *ctx)
97069a06e7bSVaclav Hapla {
97169a06e7bSVaclav Hapla   HDF5ReadCtx    h;
97269a06e7bSVaclav Hapla   PetscErrorCode ierr;
97369a06e7bSVaclav Hapla 
97469a06e7bSVaclav Hapla   PetscFunctionBegin;
97569a06e7bSVaclav Hapla   h = *ctx;
976b8ef406cSVaclav Hapla   PetscStackCallHDF5(H5Pclose,(h->plist));
97769a06e7bSVaclav Hapla   if (h->group != h->file) PetscStackCallHDF5(H5Gclose,(h->group));
97869a06e7bSVaclav Hapla   PetscStackCallHDF5(H5Sclose,(h->dataspace));
97969a06e7bSVaclav Hapla   PetscStackCallHDF5(H5Dclose,(h->dataset));
98069a06e7bSVaclav Hapla   ierr = PetscFree(*ctx);CHKERRQ(ierr);
98169a06e7bSVaclav Hapla   PetscFunctionReturn(0);
98269a06e7bSVaclav Hapla }
98369a06e7bSVaclav Hapla 
984eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadSizes_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout *map_)
98569a06e7bSVaclav Hapla {
98669a06e7bSVaclav Hapla   int            rdim, dim;
98769a06e7bSVaclav Hapla   hsize_t        dims[4];
98845e21b7fSVaclav Hapla   PetscInt       bsInd, lenInd, bs, N;
9898374c777SVaclav Hapla   PetscLayout    map;
9908374c777SVaclav Hapla   PetscErrorCode ierr;
99169a06e7bSVaclav Hapla 
99269a06e7bSVaclav Hapla   PetscFunctionBegin;
9938374c777SVaclav Hapla   if (!(*map_)) {
9948374c777SVaclav Hapla     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)viewer),map_);CHKERRQ(ierr);
9958374c777SVaclav Hapla   }
9968374c777SVaclav Hapla   map = *map_;
9979787f08bSVaclav Hapla   /* calculate expected number of dimensions */
998a5e1feadSVaclav Hapla   dim = 0;
9999568d583SVaclav Hapla   if (ctx->timestep >= 0) ++dim;
1000a5e1feadSVaclav Hapla   ++dim; /* length in blocks */
10019568d583SVaclav Hapla   if (ctx->complexVal) ++dim;
10029787f08bSVaclav Hapla   /* get actual number of dimensions in dataset */
100369a06e7bSVaclav Hapla   PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(ctx->dataspace, dims, NULL));
10049787f08bSVaclav Hapla   /* calculate expected dimension indices */
10059787f08bSVaclav Hapla   lenInd = 0;
10069568d583SVaclav Hapla   if (ctx->timestep >= 0) ++lenInd;
10079787f08bSVaclav Hapla   bsInd = lenInd + 1;
10089568d583SVaclav Hapla   ctx->dim2 = PETSC_FALSE;
10099787f08bSVaclav Hapla   if (rdim == dim) {
101045e21b7fSVaclav Hapla     bs = 1; /* support vectors stored as 1D array */
10119787f08bSVaclav Hapla   } else if (rdim == dim+1) {
101245e21b7fSVaclav Hapla     bs = (PetscInt) dims[bsInd];
10139568d583SVaclav Hapla     if (bs == 1) ctx->dim2 = PETSC_TRUE; /* vector with blocksize of 1, still stored as 2D array */
1014a5e1feadSVaclav Hapla   } else {
10159787f08bSVaclav Hapla     SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Dimension of array in file %d not %d as expected", rdim, dim);
1016a5e1feadSVaclav Hapla   }
101745e21b7fSVaclav Hapla   N = (PetscInt) dims[lenInd]*bs;
10188374c777SVaclav Hapla 
10198374c777SVaclav Hapla   /* Set Vec sizes,blocksize,and type if not already set */
10208374c777SVaclav Hapla   if (map->bs < 0) {
102145e21b7fSVaclav Hapla     ierr = PetscLayoutSetBlockSize(map, bs);CHKERRQ(ierr);
102245e21b7fSVaclav Hapla   } else if (map->bs != bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Block size of array in file is %D, not %D as expected",bs,map->bs);
10238374c777SVaclav Hapla   if (map->N < 0) {
102445e21b7fSVaclav Hapla     ierr = PetscLayoutSetSize(map, N);CHKERRQ(ierr);
102545e21b7fSVaclav Hapla   } else if (map->N != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Global size of array in file is %D, not %D as expected",N,map->N);
102669a06e7bSVaclav Hapla   PetscFunctionReturn(0);
102769a06e7bSVaclav Hapla }
102869a06e7bSVaclav Hapla 
1029eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadSelectHyperslab_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout map, hid_t *memspace)
103016f6402dSVaclav Hapla {
103116f6402dSVaclav Hapla   hsize_t        count[4], offset[4];
103216f6402dSVaclav Hapla   int            dim;
103316f6402dSVaclav Hapla   PetscInt       bs, n, low;
103416f6402dSVaclav Hapla   PetscErrorCode ierr;
103516f6402dSVaclav Hapla 
103616f6402dSVaclav Hapla   PetscFunctionBegin;
103716f6402dSVaclav Hapla   /* Compute local size and ownership range */
103816f6402dSVaclav Hapla   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
103916f6402dSVaclav Hapla   ierr = PetscLayoutGetBlockSize(map, &bs);CHKERRQ(ierr);
104016f6402dSVaclav Hapla   ierr = PetscLayoutGetLocalSize(map, &n);CHKERRQ(ierr);
104116f6402dSVaclav Hapla   ierr = PetscLayoutGetRange(map, &low, NULL);CHKERRQ(ierr);
104216f6402dSVaclav Hapla 
104316f6402dSVaclav Hapla   /* Each process defines a dataset and reads it from the hyperslab in the file */
104416f6402dSVaclav Hapla   dim  = 0;
104516f6402dSVaclav Hapla   if (ctx->timestep >= 0) {
104616f6402dSVaclav Hapla     count[dim]  = 1;
104716f6402dSVaclav Hapla     offset[dim] = ctx->timestep;
104816f6402dSVaclav Hapla     ++dim;
104916f6402dSVaclav Hapla   }
105016f6402dSVaclav Hapla   {
105116f6402dSVaclav Hapla     ierr = PetscHDF5IntCast(n/bs, &count[dim]);CHKERRQ(ierr);
105216f6402dSVaclav Hapla     ierr = PetscHDF5IntCast(low/bs, &offset[dim]);CHKERRQ(ierr);
105316f6402dSVaclav Hapla     ++dim;
105416f6402dSVaclav Hapla   }
105516f6402dSVaclav Hapla   if (bs > 1 || ctx->dim2) {
105616f6402dSVaclav Hapla     count[dim]  = bs;
105716f6402dSVaclav Hapla     offset[dim] = 0;
105816f6402dSVaclav Hapla     ++dim;
105916f6402dSVaclav Hapla   }
106016f6402dSVaclav Hapla   if (ctx->complexVal) {
106116f6402dSVaclav Hapla     count[dim]  = 2;
106216f6402dSVaclav Hapla     offset[dim] = 0;
106316f6402dSVaclav Hapla     ++dim;
106416f6402dSVaclav Hapla   }
106516f6402dSVaclav Hapla   PetscStackCallHDF5Return(*memspace,H5Screate_simple,(dim, count, NULL));
106616f6402dSVaclav Hapla   PetscStackCallHDF5(H5Sselect_hyperslab,(ctx->dataspace, H5S_SELECT_SET, offset, NULL, count, NULL));
106716f6402dSVaclav Hapla   PetscFunctionReturn(0);
106816f6402dSVaclav Hapla }
106916f6402dSVaclav Hapla 
1070eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadArray_Private(PetscViewer viewer, HDF5ReadCtx h, hid_t datatype, hid_t memspace, void *arr)
1071ef2d82ceSVaclav Hapla {
1072ef2d82ceSVaclav Hapla   PetscFunctionBegin;
1073ef2d82ceSVaclav Hapla   PetscStackCallHDF5(H5Dread,(h->dataset, datatype, memspace, h->dataspace, h->plist, arr));
1074ef2d82ceSVaclav Hapla   PetscFunctionReturn(0);
1075ef2d82ceSVaclav Hapla }
1076ef2d82ceSVaclav Hapla 
1077dad982a8SVaclav Hapla PetscErrorCode PetscViewerHDF5Load(PetscViewer viewer, const char *name, PetscLayout map, hid_t datatype, void **newarr)
1078a25c73c6SVaclav Hapla {
1079fbd37863SVaclav Hapla   HDF5ReadCtx     h=NULL;
1080fbd37863SVaclav Hapla   hid_t           memspace=0;
1081a25c73c6SVaclav Hapla   size_t          unitsize;
1082a25c73c6SVaclav Hapla   void            *arr;
1083a25c73c6SVaclav Hapla   PetscErrorCode  ierr;
1084a25c73c6SVaclav Hapla 
1085a25c73c6SVaclav Hapla   PetscFunctionBegin;
1086eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr);
1087eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr);
1088a25c73c6SVaclav Hapla   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
1089eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadSelectHyperslab_Private(viewer, h, map, &memspace);CHKERRQ(ierr);
10904fc17bcdSVaclav Hapla 
10915a89fdf4SVaclav Hapla #if defined(PETSC_USE_COMPLEX)
10925a89fdf4SVaclav Hapla   if (!h->complexVal) {
1093c76551afSVaclav Hapla     H5T_class_t clazz = H5Tget_class(datatype);
1094c76551afSVaclav Hapla     if (clazz == H5T_FLOAT) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"File contains real numbers but PETSc is configured for complex. The conversion is not yet implemented. Configure with --with-scalar-type=real.");
10955a89fdf4SVaclav Hapla   }
10965a89fdf4SVaclav Hapla #else
10975a89fdf4SVaclav Hapla   if (h->complexVal) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"File contains complex numbers but PETSc not configured for them. Configure with --with-scalar-type=complex.");
10985a89fdf4SVaclav Hapla #endif
10994fc17bcdSVaclav Hapla   unitsize = H5Tget_size(datatype);
11004fc17bcdSVaclav Hapla   if (h->complexVal) unitsize *= 2;
11014fc17bcdSVaclav Hapla   ierr = PetscMalloc(map->n*unitsize, &arr);CHKERRQ(ierr);
11024fc17bcdSVaclav Hapla 
1103eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadArray_Private(viewer, h, datatype, memspace, arr);CHKERRQ(ierr);
1104a25c73c6SVaclav Hapla   PetscStackCallHDF5(H5Sclose,(memspace));
1105eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr);
1106a25c73c6SVaclav Hapla   *newarr = arr;
1107a25c73c6SVaclav Hapla   PetscFunctionReturn(0);
1108a25c73c6SVaclav Hapla }
1109a25c73c6SVaclav Hapla 
1110c1aaad9cSVaclav Hapla /*@C
1111c1aaad9cSVaclav Hapla  PetscViewerHDF5ReadSizes - Read block size and global size of a vector (Vec or IS) stored in an HDF5 file.
1112c1aaad9cSVaclav Hapla 
1113c1aaad9cSVaclav Hapla   Input Parameters:
1114c1aaad9cSVaclav Hapla + viewer - The HDF5 viewer
1115c1aaad9cSVaclav Hapla - name   - The vector name
1116c1aaad9cSVaclav Hapla 
1117c1aaad9cSVaclav Hapla   Output Parameter:
1118c1aaad9cSVaclav Hapla + bs     - block size
1119c1aaad9cSVaclav Hapla - N      - global size
1120c1aaad9cSVaclav Hapla 
1121c1aaad9cSVaclav Hapla   Note:
1122c1aaad9cSVaclav Hapla   A vector is stored as an HDF5 dataspace with 1-4 dimensions in this order:
1123c1aaad9cSVaclav Hapla   1) # timesteps (optional), 2) # blocks, 3) # elements per block (optional), 4) real and imaginary part (only for complex).
1124c1aaad9cSVaclav Hapla 
1125c1aaad9cSVaclav Hapla   A vectors can be stored as a 2D dataspace even if its blocksize is 1; see PetscViewerHDF5SetBaseDimension2().
1126c1aaad9cSVaclav Hapla 
1127c1aaad9cSVaclav Hapla   Level: advanced
1128c1aaad9cSVaclav Hapla 
1129c1aaad9cSVaclav Hapla .seealso: PetscViewerHDF5Open(), VecLoad(), ISLoad(), VecGetSize(), ISGetSize(), PetscViewerHDF5SetBaseDimension2()
1130c1aaad9cSVaclav Hapla @*/
113169a06e7bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadSizes(PetscViewer viewer, const char name[], PetscInt *bs, PetscInt *N)
113269a06e7bSVaclav Hapla {
1133fbd37863SVaclav Hapla   HDF5ReadCtx    h=NULL;
11348374c777SVaclav Hapla   PetscLayout    map=NULL;
113569a06e7bSVaclav Hapla   PetscErrorCode ierr;
113669a06e7bSVaclav Hapla 
113769a06e7bSVaclav Hapla   PetscFunctionBegin;
1138c1aaad9cSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1139eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr);
1140eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr);
1141eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr);
11428374c777SVaclav Hapla   if (bs) *bs = map->bs;
11438374c777SVaclav Hapla   if (N) *N = map->N;
11448374c777SVaclav Hapla   ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr);
1145a5e1feadSVaclav Hapla   PetscFunctionReturn(0);
1146a5e1feadSVaclav Hapla }
1147a5e1feadSVaclav Hapla 
1148a75e6a4aSMatthew G. Knepley /*
1149a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1150a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
1151a75e6a4aSMatthew G. Knepley */
1152d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1153a75e6a4aSMatthew G. Knepley 
1154a75e6a4aSMatthew G. Knepley /*@C
1155a75e6a4aSMatthew G. Knepley   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
1156a75e6a4aSMatthew G. Knepley 
1157a75e6a4aSMatthew G. Knepley   Collective on MPI_Comm
1158a75e6a4aSMatthew G. Knepley 
1159a75e6a4aSMatthew G. Knepley   Input Parameter:
1160a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer
1161a75e6a4aSMatthew G. Knepley 
1162a75e6a4aSMatthew G. Knepley   Level: intermediate
1163a75e6a4aSMatthew G. Knepley 
1164a75e6a4aSMatthew G. Knepley   Options Database Keys:
1165a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name>
1166a75e6a4aSMatthew G. Knepley 
1167a75e6a4aSMatthew G. Knepley   Environmental variables:
1168a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME
1169a75e6a4aSMatthew G. Knepley 
1170a75e6a4aSMatthew G. Knepley   Notes:
1171a75e6a4aSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
1172a75e6a4aSMatthew G. Knepley   an error code.  The HDF5 PetscViewer is usually used in the form
1173a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1174a75e6a4aSMatthew G. Knepley 
1175a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1176a75e6a4aSMatthew G. Knepley @*/
1177a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1178a75e6a4aSMatthew G. Knepley {
1179a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
1180a75e6a4aSMatthew G. Knepley   PetscBool      flg;
1181a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
1182a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
1183a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
1184a75e6a4aSMatthew G. Knepley 
1185a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
1186a75e6a4aSMatthew 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);}
1187a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
118812801b39SBarry Smith     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
1189a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1190a75e6a4aSMatthew G. Knepley   }
119147435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1192a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1193a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
1194a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
1195a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1196a75e6a4aSMatthew G. Knepley     if (!flg) {
1197a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname,"output.h5");
1198a75e6a4aSMatthew G. Knepley       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1199a75e6a4aSMatthew G. Knepley     }
1200a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
1201a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1202a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
1203a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
120447435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
1205a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1206a75e6a4aSMatthew G. Knepley   }
1207a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
1208a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1209a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
1210a75e6a4aSMatthew G. Knepley }
1211