xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision dad982a8d69e7a6fc03bc9ef92d73188b0f31e8a)
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 
19eb5a92b4SVaclav Hapla struct _n_HDF5ReadCtx {
20eb5a92b4SVaclav Hapla   hid_t file, group, dataset, dataspace, plist;
21eb5a92b4SVaclav Hapla   PetscInt timestep;
22eb5a92b4SVaclav Hapla   PetscBool complexVal, dim2;
23eb5a92b4SVaclav Hapla };
24eb5a92b4SVaclav Hapla typedef struct _n_HDF5ReadCtx* HDF5ReadCtx;
25eb5a92b4SVaclav Hapla 
264416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v)
2782ea9b62SBarry Smith {
2882ea9b62SBarry Smith   PetscErrorCode   ierr;
2982ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;
3082ea9b62SBarry Smith 
3182ea9b62SBarry Smith   PetscFunctionBegin;
3282ea9b62SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr);
3382ea9b62SBarry Smith   ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr);
349a0502c6SHåkon Strandenes   ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr);
3582ea9b62SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3682ea9b62SBarry Smith   PetscFunctionReturn(0);
3782ea9b62SBarry Smith }
3882ea9b62SBarry Smith 
395c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
405c6c1daeSBarry Smith {
415c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
425c6c1daeSBarry Smith   PetscErrorCode   ierr;
435c6c1daeSBarry Smith 
445c6c1daeSBarry Smith   PetscFunctionBegin;
455c6c1daeSBarry Smith   ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);
46729ab6d9SBarry Smith   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
475c6c1daeSBarry Smith   PetscFunctionReturn(0);
485c6c1daeSBarry Smith }
495c6c1daeSBarry Smith 
505c6c1daeSBarry Smith PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
515c6c1daeSBarry Smith {
525c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
535c6c1daeSBarry Smith   PetscErrorCode   ierr;
545c6c1daeSBarry Smith 
555c6c1daeSBarry Smith   PetscFunctionBegin;
565c6c1daeSBarry Smith   ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr);
575c6c1daeSBarry Smith   while (hdf5->groups) {
585c6c1daeSBarry Smith     GroupList *tmp = hdf5->groups->next;
595c6c1daeSBarry Smith 
605c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups->name);CHKERRQ(ierr);
615c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups);CHKERRQ(ierr);
625c6c1daeSBarry Smith     hdf5->groups = tmp;
635c6c1daeSBarry Smith   }
645c6c1daeSBarry Smith   ierr = PetscFree(hdf5);CHKERRQ(ierr);
650b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
66d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr);
670b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
685c6c1daeSBarry Smith   PetscFunctionReturn(0);
695c6c1daeSBarry Smith }
705c6c1daeSBarry Smith 
715c6c1daeSBarry Smith PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
725c6c1daeSBarry Smith {
735c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
745c6c1daeSBarry Smith 
755c6c1daeSBarry Smith   PetscFunctionBegin;
765c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
775c6c1daeSBarry Smith   hdf5->btype = type;
785c6c1daeSBarry Smith   PetscFunctionReturn(0);
795c6c1daeSBarry Smith }
805c6c1daeSBarry Smith 
8182ea9b62SBarry Smith PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
8282ea9b62SBarry Smith {
8382ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
8482ea9b62SBarry Smith 
8582ea9b62SBarry Smith   PetscFunctionBegin;
8682ea9b62SBarry Smith   hdf5->basedimension2 = flg;
8782ea9b62SBarry Smith   PetscFunctionReturn(0);
8882ea9b62SBarry Smith }
8982ea9b62SBarry Smith 
90df863907SAlex Fikl /*@
9182ea9b62SBarry Smith      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
9282ea9b62SBarry Smith        dimension of 2.
9382ea9b62SBarry Smith 
9482ea9b62SBarry Smith     Logically Collective on PetscViewer
9582ea9b62SBarry Smith 
9682ea9b62SBarry Smith   Input Parameters:
9782ea9b62SBarry Smith +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
9882ea9b62SBarry 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
9982ea9b62SBarry Smith 
10082ea9b62SBarry Smith   Options Database:
10182ea9b62SBarry 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
10282ea9b62SBarry Smith 
10382ea9b62SBarry Smith 
10495452b02SPatrick Sanan   Notes:
10595452b02SPatrick 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
10682ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
10782ea9b62SBarry Smith 
10882ea9b62SBarry Smith   Level: intermediate
10982ea9b62SBarry Smith 
11082ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
11182ea9b62SBarry Smith 
11282ea9b62SBarry Smith @*/
11382ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
11482ea9b62SBarry Smith {
11582ea9b62SBarry Smith   PetscErrorCode ierr;
11682ea9b62SBarry Smith 
11782ea9b62SBarry Smith   PetscFunctionBegin;
11882ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
11982ea9b62SBarry Smith   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
12082ea9b62SBarry Smith   PetscFunctionReturn(0);
12182ea9b62SBarry Smith }
12282ea9b62SBarry Smith 
123df863907SAlex Fikl /*@
12482ea9b62SBarry Smith      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
12582ea9b62SBarry Smith        dimension of 2.
12682ea9b62SBarry Smith 
12782ea9b62SBarry Smith     Logically Collective on PetscViewer
12882ea9b62SBarry Smith 
12982ea9b62SBarry Smith   Input Parameter:
13082ea9b62SBarry Smith .  viewer - the PetscViewer, must be of type HDF5
13182ea9b62SBarry Smith 
13282ea9b62SBarry Smith   Output Parameter:
13382ea9b62SBarry 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
13482ea9b62SBarry Smith 
13595452b02SPatrick Sanan   Notes:
13695452b02SPatrick 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
13782ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
13882ea9b62SBarry Smith 
13982ea9b62SBarry Smith   Level: intermediate
14082ea9b62SBarry Smith 
14182ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
14282ea9b62SBarry Smith 
14382ea9b62SBarry Smith @*/
14482ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
14582ea9b62SBarry Smith {
14682ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
14782ea9b62SBarry Smith 
14882ea9b62SBarry Smith   PetscFunctionBegin;
14982ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
15082ea9b62SBarry Smith   *flg = hdf5->basedimension2;
15182ea9b62SBarry Smith   PetscFunctionReturn(0);
15282ea9b62SBarry Smith }
15382ea9b62SBarry Smith 
1549a0502c6SHåkon Strandenes PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
1559a0502c6SHåkon Strandenes {
1569a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1579a0502c6SHåkon Strandenes 
1589a0502c6SHåkon Strandenes   PetscFunctionBegin;
1599a0502c6SHåkon Strandenes   hdf5->spoutput = flg;
1609a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
1619a0502c6SHåkon Strandenes }
1629a0502c6SHåkon Strandenes 
163df863907SAlex Fikl /*@
1649a0502c6SHåkon Strandenes      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
1659a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
1669a0502c6SHåkon Strandenes 
1679a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
1689a0502c6SHåkon Strandenes 
1699a0502c6SHåkon Strandenes   Input Parameters:
1709a0502c6SHåkon Strandenes +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
1719a0502c6SHåkon Strandenes -  flg - if PETSC_TRUE the data will be written to disk with single precision
1729a0502c6SHåkon Strandenes 
1739a0502c6SHåkon Strandenes   Options Database:
1749a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
1759a0502c6SHåkon Strandenes 
1769a0502c6SHåkon Strandenes 
17795452b02SPatrick Sanan   Notes:
17895452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
1799a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
1809a0502c6SHåkon Strandenes 
1819a0502c6SHåkon Strandenes   Level: intermediate
1829a0502c6SHåkon Strandenes 
1839a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
1849a0502c6SHåkon Strandenes           PetscReal
1859a0502c6SHåkon Strandenes 
1869a0502c6SHåkon Strandenes @*/
1879a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
1889a0502c6SHåkon Strandenes {
1899a0502c6SHåkon Strandenes   PetscErrorCode ierr;
1909a0502c6SHåkon Strandenes 
1919a0502c6SHåkon Strandenes   PetscFunctionBegin;
1929a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1939a0502c6SHåkon Strandenes   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
1949a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
1959a0502c6SHåkon Strandenes }
1969a0502c6SHåkon Strandenes 
197df863907SAlex Fikl /*@
1989a0502c6SHåkon Strandenes      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
1999a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2009a0502c6SHåkon Strandenes 
2019a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2029a0502c6SHåkon Strandenes 
2039a0502c6SHåkon Strandenes   Input Parameter:
2049a0502c6SHåkon Strandenes .  viewer - the PetscViewer, must be of type HDF5
2059a0502c6SHåkon Strandenes 
2069a0502c6SHåkon Strandenes   Output Parameter:
2079a0502c6SHåkon Strandenes .  flg - if PETSC_TRUE the data will be written to disk with single precision
2089a0502c6SHåkon Strandenes 
20995452b02SPatrick Sanan   Notes:
21095452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2119a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2129a0502c6SHåkon Strandenes 
2139a0502c6SHåkon Strandenes   Level: intermediate
2149a0502c6SHåkon Strandenes 
2159a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2169a0502c6SHåkon Strandenes           PetscReal
2179a0502c6SHåkon Strandenes 
2189a0502c6SHåkon Strandenes @*/
2199a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
2209a0502c6SHåkon Strandenes {
2219a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2229a0502c6SHåkon Strandenes 
2239a0502c6SHåkon Strandenes   PetscFunctionBegin;
2249a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2259a0502c6SHåkon Strandenes   *flg = hdf5->spoutput;
2269a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2279a0502c6SHåkon Strandenes }
2289a0502c6SHåkon Strandenes 
2295c6c1daeSBarry Smith PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
2305c6c1daeSBarry Smith {
2315c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2325c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
2335c6c1daeSBarry Smith   MPI_Info          info = MPI_INFO_NULL;
2345c6c1daeSBarry Smith #endif
2355c6c1daeSBarry Smith   hid_t             plist_id;
2365c6c1daeSBarry Smith   PetscErrorCode    ierr;
2375c6c1daeSBarry Smith 
2385c6c1daeSBarry Smith   PetscFunctionBegin;
239f683cc58SToby Isaac   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
240f683cc58SToby Isaac   if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);}
2415c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
2425c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
243729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
2445c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
245729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info));
2465c6c1daeSBarry Smith #endif
2475c6c1daeSBarry Smith   /* Create or open the file collectively */
2485c6c1daeSBarry Smith   switch (hdf5->btype) {
2495c6c1daeSBarry Smith   case FILE_MODE_READ:
250729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
2515c6c1daeSBarry Smith     break;
2525c6c1daeSBarry Smith   case FILE_MODE_APPEND:
253729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
2545c6c1daeSBarry Smith     break;
2555c6c1daeSBarry Smith   case FILE_MODE_WRITE:
256729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
2575c6c1daeSBarry Smith     break;
2585c6c1daeSBarry Smith   default:
2595c6c1daeSBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
2605c6c1daeSBarry Smith   }
2615c6c1daeSBarry Smith   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
262729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
2635c6c1daeSBarry Smith   PetscFunctionReturn(0);
2645c6c1daeSBarry Smith }
2655c6c1daeSBarry Smith 
266d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
267d1232d7fSToby Isaac {
268d1232d7fSToby Isaac   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
269d1232d7fSToby Isaac 
270d1232d7fSToby Isaac   PetscFunctionBegin;
271d1232d7fSToby Isaac   *name = vhdf5->filename;
272d1232d7fSToby Isaac   PetscFunctionReturn(0);
273d1232d7fSToby Isaac }
274d1232d7fSToby Isaac 
2758556b5ebSBarry Smith /*MC
2768556b5ebSBarry Smith    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
2778556b5ebSBarry Smith 
2788556b5ebSBarry Smith 
2798556b5ebSBarry Smith .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
2808556b5ebSBarry Smith            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
2818556b5ebSBarry Smith            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
2828556b5ebSBarry Smith            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
2838556b5ebSBarry Smith 
2841b266c99SBarry Smith   Level: beginner
2858556b5ebSBarry Smith M*/
286d1232d7fSToby Isaac 
2878cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
2885c6c1daeSBarry Smith {
2895c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
2905c6c1daeSBarry Smith   PetscErrorCode   ierr;
2915c6c1daeSBarry Smith 
2925c6c1daeSBarry Smith   PetscFunctionBegin;
293b00a9115SJed Brown   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
2945c6c1daeSBarry Smith 
2955c6c1daeSBarry Smith   v->data                = (void*) hdf5;
2965c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
29782ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
2985c6c1daeSBarry Smith   v->ops->flush          = 0;
2995c6c1daeSBarry Smith   hdf5->btype            = (PetscFileMode) -1;
3005c6c1daeSBarry Smith   hdf5->filename         = 0;
3015c6c1daeSBarry Smith   hdf5->timestep         = -1;
3020298fd71SBarry Smith   hdf5->groups           = NULL;
3035c6c1daeSBarry Smith 
3040b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
305d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr);
3060b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
30782ea9b62SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
3089a0502c6SHåkon Strandenes   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr);
3095c6c1daeSBarry Smith   PetscFunctionReturn(0);
3105c6c1daeSBarry Smith }
3115c6c1daeSBarry Smith 
3125c6c1daeSBarry Smith /*@C
3135c6c1daeSBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
3145c6c1daeSBarry Smith 
3155c6c1daeSBarry Smith    Collective on MPI_Comm
3165c6c1daeSBarry Smith 
3175c6c1daeSBarry Smith    Input Parameters:
3185c6c1daeSBarry Smith +  comm - MPI communicator
3195c6c1daeSBarry Smith .  name - name of file
3205c6c1daeSBarry Smith -  type - type of file
3215c6c1daeSBarry Smith $    FILE_MODE_WRITE - create new file for binary output
3225c6c1daeSBarry Smith $    FILE_MODE_READ - open existing file for binary input
3235c6c1daeSBarry Smith $    FILE_MODE_APPEND - open existing file for binary output
3245c6c1daeSBarry Smith 
3255c6c1daeSBarry Smith    Output Parameter:
3265c6c1daeSBarry Smith .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
3275c6c1daeSBarry Smith 
32882ea9b62SBarry Smith   Options Database:
32982ea9b62SBarry 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
3309a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
33182ea9b62SBarry Smith 
3325c6c1daeSBarry Smith    Level: beginner
3335c6c1daeSBarry Smith 
3345c6c1daeSBarry Smith    Note:
3355c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
3365c6c1daeSBarry Smith 
3375c6c1daeSBarry Smith    Concepts: HDF5 files
3385c6c1daeSBarry Smith    Concepts: PetscViewerHDF5^creating
3395c6c1daeSBarry Smith 
3406a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
3419a0502c6SHåkon Strandenes           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
342a56f64adSBarry Smith           MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
3435c6c1daeSBarry Smith @*/
3445c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
3455c6c1daeSBarry Smith {
3465c6c1daeSBarry Smith   PetscErrorCode ierr;
3475c6c1daeSBarry Smith 
3485c6c1daeSBarry Smith   PetscFunctionBegin;
3495c6c1daeSBarry Smith   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
3505c6c1daeSBarry Smith   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
3515c6c1daeSBarry Smith   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
3525c6c1daeSBarry Smith   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
3535c6c1daeSBarry Smith   PetscFunctionReturn(0);
3545c6c1daeSBarry Smith }
3555c6c1daeSBarry Smith 
3565c6c1daeSBarry Smith /*@C
3575c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
3585c6c1daeSBarry Smith 
3595c6c1daeSBarry Smith   Not collective
3605c6c1daeSBarry Smith 
3615c6c1daeSBarry Smith   Input Parameter:
3625c6c1daeSBarry Smith . viewer - the PetscViewer
3635c6c1daeSBarry Smith 
3645c6c1daeSBarry Smith   Output Parameter:
3655c6c1daeSBarry Smith . file_id - The file id
3665c6c1daeSBarry Smith 
3675c6c1daeSBarry Smith   Level: intermediate
3685c6c1daeSBarry Smith 
3695c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open()
3705c6c1daeSBarry Smith @*/
3715c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
3725c6c1daeSBarry Smith {
3735c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3745c6c1daeSBarry Smith 
3755c6c1daeSBarry Smith   PetscFunctionBegin;
3765c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
3775c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
3785c6c1daeSBarry Smith   PetscFunctionReturn(0);
3795c6c1daeSBarry Smith }
3805c6c1daeSBarry Smith 
3815c6c1daeSBarry Smith /*@C
3825c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
3835c6c1daeSBarry Smith 
3845c6c1daeSBarry Smith   Not collective
3855c6c1daeSBarry Smith 
3865c6c1daeSBarry Smith   Input Parameters:
3875c6c1daeSBarry Smith + viewer - the PetscViewer
3885c6c1daeSBarry Smith - name - The group name
3895c6c1daeSBarry Smith 
3905c6c1daeSBarry Smith   Level: intermediate
3915c6c1daeSBarry Smith 
392874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
3935c6c1daeSBarry Smith @*/
3945c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
3955c6c1daeSBarry Smith {
3965c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3975c6c1daeSBarry Smith   GroupList        *groupNode;
3985c6c1daeSBarry Smith   PetscErrorCode   ierr;
3995c6c1daeSBarry Smith 
4005c6c1daeSBarry Smith   PetscFunctionBegin;
4015c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
4025c6c1daeSBarry Smith   PetscValidCharPointer(name,2);
40395dccacaSBarry Smith   ierr = PetscNew(&groupNode);CHKERRQ(ierr);
4045c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
405a297a907SKarl Rupp 
4065c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
4075c6c1daeSBarry Smith   hdf5->groups    = groupNode;
4085c6c1daeSBarry Smith   PetscFunctionReturn(0);
4095c6c1daeSBarry Smith }
4105c6c1daeSBarry Smith 
4113ef9c667SSatish Balay /*@
4125c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
4135c6c1daeSBarry Smith 
4145c6c1daeSBarry Smith   Not collective
4155c6c1daeSBarry Smith 
4165c6c1daeSBarry Smith   Input Parameter:
4175c6c1daeSBarry Smith . viewer - the PetscViewer
4185c6c1daeSBarry Smith 
4195c6c1daeSBarry Smith   Level: intermediate
4205c6c1daeSBarry Smith 
421874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
4225c6c1daeSBarry Smith @*/
4235c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
4245c6c1daeSBarry Smith {
4255c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4265c6c1daeSBarry Smith   GroupList        *groupNode;
4275c6c1daeSBarry Smith   PetscErrorCode   ierr;
4285c6c1daeSBarry Smith 
4295c6c1daeSBarry Smith   PetscFunctionBegin;
4305c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
43182f516ccSBarry Smith   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
4325c6c1daeSBarry Smith   groupNode    = hdf5->groups;
4335c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
4345c6c1daeSBarry Smith   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
4355c6c1daeSBarry Smith   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
4365c6c1daeSBarry Smith   PetscFunctionReturn(0);
4375c6c1daeSBarry Smith }
4385c6c1daeSBarry Smith 
4395c6c1daeSBarry Smith /*@C
440874270d9SVaclav Hapla   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup().
441874270d9SVaclav Hapla   If none has been assigned, returns NULL.
4425c6c1daeSBarry Smith 
4435c6c1daeSBarry Smith   Not collective
4445c6c1daeSBarry Smith 
4455c6c1daeSBarry Smith   Input Parameter:
4465c6c1daeSBarry Smith . viewer - the PetscViewer
4475c6c1daeSBarry Smith 
4485c6c1daeSBarry Smith   Output Parameter:
4495c6c1daeSBarry Smith . name - The group name
4505c6c1daeSBarry Smith 
4515c6c1daeSBarry Smith   Level: intermediate
4525c6c1daeSBarry Smith 
453874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
4545c6c1daeSBarry Smith @*/
4555c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
4565c6c1daeSBarry Smith {
4575c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
4585c6c1daeSBarry Smith 
4595c6c1daeSBarry Smith   PetscFunctionBegin;
4605c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
461c959eef4SJed Brown   PetscValidPointer(name,2);
462a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
4630298fd71SBarry Smith   else *name = NULL;
4645c6c1daeSBarry Smith   PetscFunctionReturn(0);
4655c6c1daeSBarry Smith }
4665c6c1daeSBarry Smith 
4678c557b5aSVaclav Hapla /*@
468874270d9SVaclav Hapla   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(),
469874270d9SVaclav Hapla   and return this group's ID and file ID.
470874270d9SVaclav Hapla   If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID.
471874270d9SVaclav Hapla 
472874270d9SVaclav Hapla   Not collective
473874270d9SVaclav Hapla 
474874270d9SVaclav Hapla   Input Parameter:
475874270d9SVaclav Hapla . viewer - the PetscViewer
476874270d9SVaclav Hapla 
477874270d9SVaclav Hapla   Output Parameter:
478874270d9SVaclav Hapla + fileId - The HDF5 file ID
479874270d9SVaclav Hapla - groupId - The HDF5 group ID
480874270d9SVaclav Hapla 
481874270d9SVaclav Hapla   Level: intermediate
482874270d9SVaclav Hapla 
483874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
484874270d9SVaclav Hapla @*/
48554dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
48654dbf706SVaclav Hapla {
48754dbf706SVaclav Hapla   hid_t          file_id, group;
48854dbf706SVaclav Hapla   htri_t         found;
48954dbf706SVaclav Hapla   const char     *groupName = NULL;
49054dbf706SVaclav Hapla   PetscErrorCode ierr;
49154dbf706SVaclav Hapla 
49254dbf706SVaclav Hapla   PetscFunctionBegin;
49354dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr);
49454dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr);
49554dbf706SVaclav Hapla   /* Open group */
49654dbf706SVaclav Hapla   if (groupName) {
49754dbf706SVaclav Hapla     PetscBool root;
49854dbf706SVaclav Hapla 
49954dbf706SVaclav Hapla     ierr = PetscStrcmp(groupName, "/", &root);CHKERRQ(ierr);
50054dbf706SVaclav Hapla     PetscStackCall("H5Lexists",found = H5Lexists(file_id, groupName, H5P_DEFAULT));
50154dbf706SVaclav Hapla     if (!root && (found <= 0)) {
50254dbf706SVaclav Hapla #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
50354dbf706SVaclav Hapla       PetscStackCallHDF5Return(group,H5Gcreate2,(file_id, groupName, 0, H5P_DEFAULT, H5P_DEFAULT));
50454dbf706SVaclav Hapla #else /* deprecated HDF5 1.6 API */
50554dbf706SVaclav Hapla       PetscStackCallHDF5Return(group,H5Gcreate,(file_id, groupName, 0));
50654dbf706SVaclav Hapla #endif
50754dbf706SVaclav Hapla       PetscStackCallHDF5(H5Gclose,(group));
50854dbf706SVaclav Hapla     }
50954dbf706SVaclav Hapla #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
51054dbf706SVaclav Hapla     PetscStackCallHDF5Return(group,H5Gopen2,(file_id, groupName, H5P_DEFAULT));
51154dbf706SVaclav Hapla #else
51254dbf706SVaclav Hapla     PetscStackCallHDF5Return(group,H5Gopen,(file_id, groupName));
51354dbf706SVaclav Hapla #endif
51454dbf706SVaclav Hapla   } else group = file_id;
51554dbf706SVaclav Hapla 
51654dbf706SVaclav Hapla   *fileId  = file_id;
51754dbf706SVaclav Hapla   *groupId = group;
51854dbf706SVaclav Hapla   PetscFunctionReturn(0);
51954dbf706SVaclav Hapla }
52054dbf706SVaclav Hapla 
5213ef9c667SSatish Balay /*@
5225c6c1daeSBarry Smith   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
5235c6c1daeSBarry Smith 
5245c6c1daeSBarry Smith   Not collective
5255c6c1daeSBarry Smith 
5265c6c1daeSBarry Smith   Input Parameter:
5275c6c1daeSBarry Smith . viewer - the PetscViewer
5285c6c1daeSBarry Smith 
5295c6c1daeSBarry Smith   Level: intermediate
5305c6c1daeSBarry Smith 
5315c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
5325c6c1daeSBarry Smith @*/
5335c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
5345c6c1daeSBarry Smith {
5355c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5365c6c1daeSBarry Smith 
5375c6c1daeSBarry Smith   PetscFunctionBegin;
5385c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5395c6c1daeSBarry Smith   ++hdf5->timestep;
5405c6c1daeSBarry Smith   PetscFunctionReturn(0);
5415c6c1daeSBarry Smith }
5425c6c1daeSBarry Smith 
5433ef9c667SSatish Balay /*@
5445c6c1daeSBarry Smith   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
5455c6c1daeSBarry Smith   of -1 disables blocking with timesteps.
5465c6c1daeSBarry Smith 
5475c6c1daeSBarry Smith   Not collective
5485c6c1daeSBarry Smith 
5495c6c1daeSBarry Smith   Input Parameters:
5505c6c1daeSBarry Smith + viewer - the PetscViewer
5515c6c1daeSBarry Smith - timestep - The timestep number
5525c6c1daeSBarry Smith 
5535c6c1daeSBarry Smith   Level: intermediate
5545c6c1daeSBarry Smith 
5555c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
5565c6c1daeSBarry Smith @*/
5575c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
5585c6c1daeSBarry Smith {
5595c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5605c6c1daeSBarry Smith 
5615c6c1daeSBarry Smith   PetscFunctionBegin;
5625c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5635c6c1daeSBarry Smith   hdf5->timestep = timestep;
5645c6c1daeSBarry Smith   PetscFunctionReturn(0);
5655c6c1daeSBarry Smith }
5665c6c1daeSBarry Smith 
5673ef9c667SSatish Balay /*@
5685c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
5695c6c1daeSBarry Smith 
5705c6c1daeSBarry Smith   Not collective
5715c6c1daeSBarry Smith 
5725c6c1daeSBarry Smith   Input Parameter:
5735c6c1daeSBarry Smith . viewer - the PetscViewer
5745c6c1daeSBarry Smith 
5755c6c1daeSBarry Smith   Output Parameter:
5765c6c1daeSBarry Smith . timestep - The timestep number
5775c6c1daeSBarry Smith 
5785c6c1daeSBarry Smith   Level: intermediate
5795c6c1daeSBarry Smith 
5805c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
5815c6c1daeSBarry Smith @*/
5825c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
5835c6c1daeSBarry Smith {
5845c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5855c6c1daeSBarry Smith 
5865c6c1daeSBarry Smith   PetscFunctionBegin;
5875c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5885c6c1daeSBarry Smith   PetscValidPointer(timestep,2);
5895c6c1daeSBarry Smith   *timestep = hdf5->timestep;
5905c6c1daeSBarry Smith   PetscFunctionReturn(0);
5915c6c1daeSBarry Smith }
5925c6c1daeSBarry Smith 
59336bce6e8SMatthew G. Knepley /*@C
59436bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
59536bce6e8SMatthew G. Knepley 
59636bce6e8SMatthew G. Knepley   Not collective
59736bce6e8SMatthew G. Knepley 
59836bce6e8SMatthew G. Knepley   Input Parameter:
59936bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
60036bce6e8SMatthew G. Knepley 
60136bce6e8SMatthew G. Knepley   Output Parameter:
60236bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
60336bce6e8SMatthew G. Knepley 
60436bce6e8SMatthew G. Knepley   Level: advanced
60536bce6e8SMatthew G. Knepley 
60636bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
60736bce6e8SMatthew G. Knepley @*/
60836bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
60936bce6e8SMatthew G. Knepley {
61036bce6e8SMatthew G. Knepley   PetscFunctionBegin;
61136bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
61236bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
61336bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_LLONG;
61436bce6e8SMatthew G. Knepley #else
61536bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_INT;
61636bce6e8SMatthew G. Knepley #endif
61736bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
61836bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
61936bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
62036bce6e8SMatthew G. Knepley   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_DOUBLE;
62136bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_DOUBLE;
62236bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
62336bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
62436bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
6257e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
62636bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
62736bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
62836bce6e8SMatthew G. Knepley }
62936bce6e8SMatthew G. Knepley 
63036bce6e8SMatthew G. Knepley /*@C
63136bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
63236bce6e8SMatthew G. Knepley 
63336bce6e8SMatthew G. Knepley   Not collective
63436bce6e8SMatthew G. Knepley 
63536bce6e8SMatthew G. Knepley   Input Parameter:
63636bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
63736bce6e8SMatthew G. Knepley 
63836bce6e8SMatthew G. Knepley   Output Parameter:
63936bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
64036bce6e8SMatthew G. Knepley 
64136bce6e8SMatthew G. Knepley   Level: advanced
64236bce6e8SMatthew G. Knepley 
64336bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
64436bce6e8SMatthew G. Knepley @*/
64536bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
64636bce6e8SMatthew G. Knepley {
64736bce6e8SMatthew G. Knepley   PetscFunctionBegin;
64836bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
64936bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
65036bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
65136bce6e8SMatthew G. Knepley #else
65236bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
65336bce6e8SMatthew G. Knepley #endif
65436bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
65536bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
65636bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
65736bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
65836bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
65936bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
6607e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
66136bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
66236bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
66336bce6e8SMatthew G. Knepley }
66436bce6e8SMatthew G. Knepley 
665df863907SAlex Fikl /*@C
66636bce6e8SMatthew G. Knepley  PetscViewerHDF5WriteAttribute - Write a scalar attribute
66736bce6e8SMatthew G. Knepley 
66836bce6e8SMatthew G. Knepley   Input Parameters:
66936bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer
67036bce6e8SMatthew G. Knepley . parent - The parent name
67136bce6e8SMatthew G. Knepley . name   - The attribute name
67236bce6e8SMatthew G. Knepley . datatype - The attribute type
67336bce6e8SMatthew G. Knepley - value    - The attribute value
67436bce6e8SMatthew G. Knepley 
67536bce6e8SMatthew G. Knepley   Level: advanced
67636bce6e8SMatthew G. Knepley 
677e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute()
67836bce6e8SMatthew G. Knepley @*/
67936bce6e8SMatthew G. Knepley PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
68036bce6e8SMatthew G. Knepley {
68160568592SMatthew G. Knepley   hid_t          h5, dataspace, obj, attribute, dtype;
68236bce6e8SMatthew G. Knepley   PetscErrorCode ierr;
68336bce6e8SMatthew G. Knepley 
68436bce6e8SMatthew G. Knepley   PetscFunctionBegin;
6855cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
68636bce6e8SMatthew G. Knepley   PetscValidPointer(parent, 2);
68736bce6e8SMatthew G. Knepley   PetscValidPointer(name, 3);
68836bce6e8SMatthew G. Knepley   PetscValidPointer(value, 4);
68936bce6e8SMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
6907e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
6917e97c476SMatthew G. Knepley     size_t len;
6927e97c476SMatthew G. Knepley     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
693729ab6d9SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
6947e97c476SMatthew G. Knepley   }
69536bce6e8SMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
696729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
69760568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
69836bce6e8SMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
69960568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
70036bce6e8SMatthew G. Knepley #else
70160568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Acreate,(obj, name, dtype, dataspace, H5P_DEFAULT));
70236bce6e8SMatthew G. Knepley #endif
703729ab6d9SBarry Smith   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
704729ab6d9SBarry Smith   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
705729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
70660568592SMatthew G. Knepley   PetscStackCallHDF5(H5Oclose,(obj));
707729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(dataspace));
70836bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
70936bce6e8SMatthew G. Knepley }
71036bce6e8SMatthew G. Knepley 
711df863907SAlex Fikl /*@C
712ad0c61feSMatthew G. Knepley  PetscViewerHDF5ReadAttribute - Read a scalar attribute
713ad0c61feSMatthew G. Knepley 
714ad0c61feSMatthew G. Knepley   Input Parameters:
715ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer
716ad0c61feSMatthew G. Knepley . parent - The parent name
717ad0c61feSMatthew G. Knepley . name   - The attribute name
718ad0c61feSMatthew G. Knepley - datatype - The attribute type
719ad0c61feSMatthew G. Knepley 
720ad0c61feSMatthew G. Knepley   Output Parameter:
721ad0c61feSMatthew G. Knepley . value    - The attribute value
722ad0c61feSMatthew G. Knepley 
723ad0c61feSMatthew G. Knepley   Level: advanced
724ad0c61feSMatthew G. Knepley 
725e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute()
726ad0c61feSMatthew G. Knepley @*/
727ad0c61feSMatthew G. Knepley PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, void *value)
728ad0c61feSMatthew G. Knepley {
729f0b58479SMatthew G. Knepley   hid_t          h5, obj, attribute, atype, dtype;
730ad0c61feSMatthew G. Knepley   PetscErrorCode ierr;
731ad0c61feSMatthew G. Knepley 
732ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
7335cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
734ad0c61feSMatthew G. Knepley   PetscValidPointer(parent, 2);
735ad0c61feSMatthew G. Knepley   PetscValidPointer(name, 3);
736ad0c61feSMatthew G. Knepley   PetscValidPointer(value, 4);
737ad0c61feSMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
738ad0c61feSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
73960568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
74060568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
741f0b58479SMatthew G. Knepley   PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
742f0b58479SMatthew G. Knepley   if (datatype == PETSC_STRING) {
743f0b58479SMatthew G. Knepley     size_t len;
744f0b58479SMatthew G. Knepley 
745f0b58479SMatthew G. Knepley     PetscStackCallHDF5Return(len,H5Tget_size,(atype));
746f0b58479SMatthew G. Knepley     PetscStackCallHDF5(H5Tclose,(atype));
747f0b58479SMatthew G. Knepley     ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr);
748f0b58479SMatthew G. Knepley   }
74970efba33SBarry Smith   PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
750729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
75160568592SMatthew G. Knepley   PetscStackCallHDF5(H5Dclose,(obj));
752ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
753ad0c61feSMatthew G. Knepley }
754ad0c61feSMatthew G. Knepley 
7555cdeb108SMatthew G. Knepley static PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, const char name[], H5O_type_t otype, PetscBool *has)
7565cdeb108SMatthew G. Knepley {
7575cdeb108SMatthew G. Knepley   hid_t          h5;
7585cdeb108SMatthew G. Knepley   PetscErrorCode ierr;
7595cdeb108SMatthew G. Knepley 
7605cdeb108SMatthew G. Knepley   PetscFunctionBegin;
7615cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7625cdeb108SMatthew G. Knepley   PetscValidPointer(name, 2);
7635cdeb108SMatthew G. Knepley   PetscValidPointer(has, 3);
7645cdeb108SMatthew G. Knepley   *has = PETSC_FALSE;
765c0cd0301SJed Brown   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
7665cdeb108SMatthew G. Knepley   if (H5Lexists(h5, name, H5P_DEFAULT)) {
7675cdeb108SMatthew G. Knepley     H5O_info_t info;
7685cdeb108SMatthew G. Knepley     hid_t      obj;
7695cdeb108SMatthew G. Knepley 
770729ab6d9SBarry Smith     PetscStackCallHDF5Return(obj,H5Oopen,(h5, name, H5P_DEFAULT));
771729ab6d9SBarry Smith     PetscStackCallHDF5(H5Oget_info,(obj, &info));
7725cdeb108SMatthew G. Knepley     if (otype == info.type) *has = PETSC_TRUE;
773729ab6d9SBarry Smith     PetscStackCallHDF5(H5Oclose,(obj));
7745cdeb108SMatthew G. Knepley   }
7755cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
7765cdeb108SMatthew G. Knepley }
7775cdeb108SMatthew G. Knepley 
778df863907SAlex Fikl /*@C
779e7bdbf8eSMatthew G. Knepley  PetscViewerHDF5HasAttribute - Check whether a scalar attribute exists
780e7bdbf8eSMatthew G. Knepley 
781e7bdbf8eSMatthew G. Knepley   Input Parameters:
782e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer
783e7bdbf8eSMatthew G. Knepley . parent - The parent name
784e7bdbf8eSMatthew G. Knepley - name   - The attribute name
785e7bdbf8eSMatthew G. Knepley 
786e7bdbf8eSMatthew G. Knepley   Output Parameter:
787e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
788e7bdbf8eSMatthew G. Knepley 
789e7bdbf8eSMatthew G. Knepley   Level: advanced
790e7bdbf8eSMatthew G. Knepley 
791e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute()
792e7bdbf8eSMatthew G. Knepley @*/
793e7bdbf8eSMatthew G. Knepley PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
794e7bdbf8eSMatthew G. Knepley {
795729ab6d9SBarry Smith   hid_t          h5, dataset;
7965cdeb108SMatthew G. Knepley   htri_t         hhas;
7975cdeb108SMatthew G. Knepley   PetscBool      exists;
798e7bdbf8eSMatthew G. Knepley   PetscErrorCode ierr;
799e7bdbf8eSMatthew G. Knepley 
800e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
8015cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
802e7bdbf8eSMatthew G. Knepley   PetscValidPointer(parent, 2);
803e7bdbf8eSMatthew G. Knepley   PetscValidPointer(name, 3);
804e7bdbf8eSMatthew G. Knepley   PetscValidPointer(has, 4);
805e7bdbf8eSMatthew G. Knepley   *has = PETSC_FALSE;
806e7bdbf8eSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
8075cdeb108SMatthew G. Knepley   ierr = PetscViewerHDF5HasObject(viewer, parent, H5O_TYPE_DATASET, &exists);CHKERRQ(ierr);
8085cdeb108SMatthew G. Knepley   if (exists) {
809e7bdbf8eSMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
810729ab6d9SBarry Smith     PetscStackCall("H5Dopen2",dataset = H5Dopen2(h5, parent, H5P_DEFAULT));
811e7bdbf8eSMatthew G. Knepley #else
812729ab6d9SBarry Smith     PetscStackCall("H5Dopen",dataset = H5Dopen(h5, parent));
813e7bdbf8eSMatthew G. Knepley #endif
814729ab6d9SBarry Smith     if (dataset < 0) PetscFunctionReturn(0);
815729ab6d9SBarry Smith     PetscStackCall("H5Aexists",hhas = H5Aexists(dataset, name));
816729ab6d9SBarry Smith     if (hhas < 0) {
817729ab6d9SBarry Smith       PetscStackCallHDF5(H5Dclose,(dataset));
818729ab6d9SBarry Smith       PetscFunctionReturn(0);
819729ab6d9SBarry Smith     }
820729ab6d9SBarry Smith     PetscStackCallHDF5(H5Dclose,(dataset));
8215cdeb108SMatthew G. Knepley     *has = hhas ? PETSC_TRUE : PETSC_FALSE;
8225cdeb108SMatthew G. Knepley   }
823e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
824e7bdbf8eSMatthew G. Knepley }
825e7bdbf8eSMatthew G. Knepley 
826eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadInitialize_Private(PetscViewer viewer, const char name[], HDF5ReadCtx *ctx)
827a5e1feadSVaclav Hapla {
82869a06e7bSVaclav Hapla   HDF5ReadCtx    h;
82969a06e7bSVaclav Hapla   const char    *groupname;
83069a06e7bSVaclav Hapla   char           vecgroup[PETSC_MAX_PATH_LEN];
831a5e1feadSVaclav Hapla   PetscErrorCode ierr;
832a5e1feadSVaclav Hapla 
833a5e1feadSVaclav Hapla   PetscFunctionBegin;
83469a06e7bSVaclav Hapla   ierr = PetscNew(&h);CHKERRQ(ierr);
83569a06e7bSVaclav Hapla   ierr = PetscViewerHDF5OpenGroup(viewer, &h->file, &h->group);CHKERRQ(ierr);
836a5e1feadSVaclav Hapla #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
83769a06e7bSVaclav Hapla   PetscStackCallHDF5Return(h->dataset,H5Dopen2,(h->group, name, H5P_DEFAULT));
838a5e1feadSVaclav Hapla #else
83969a06e7bSVaclav Hapla   PetscStackCallHDF5Return(h->dataset,H5Dopen,(h->group, name));
840a5e1feadSVaclav Hapla #endif
84169a06e7bSVaclav Hapla   PetscStackCallHDF5Return(h->dataspace,H5Dget_space,(h->dataset));
8429568d583SVaclav Hapla   ierr = PetscViewerHDF5GetTimestep(viewer, &h->timestep);CHKERRQ(ierr);
84369a06e7bSVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer,&groupname);CHKERRQ(ierr);
844b44ade6dSVaclav Hapla   ierr = PetscSNPrintf(vecgroup,PETSC_MAX_PATH_LEN,"%s/%s",groupname ? groupname : "",name);CHKERRQ(ierr);
8459568d583SVaclav Hapla   ierr = PetscViewerHDF5HasAttribute(viewer,vecgroup,"complex",&h->complexVal);CHKERRQ(ierr);
846b8ef406cSVaclav Hapla   /* Create property list for collective dataset read */
847b8ef406cSVaclav Hapla   PetscStackCallHDF5Return(h->plist,H5Pcreate,(H5P_DATASET_XFER));
848b8ef406cSVaclav Hapla #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
849b8ef406cSVaclav Hapla   PetscStackCallHDF5(H5Pset_dxpl_mpio,(h->plist, H5FD_MPIO_COLLECTIVE));
850b8ef406cSVaclav Hapla #endif
851b8ef406cSVaclav Hapla   /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */
85269a06e7bSVaclav Hapla   *ctx = h;
85369a06e7bSVaclav Hapla   PetscFunctionReturn(0);
85469a06e7bSVaclav Hapla }
85569a06e7bSVaclav Hapla 
856eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadFinalize_Private(PetscViewer viewer, HDF5ReadCtx *ctx)
85769a06e7bSVaclav Hapla {
85869a06e7bSVaclav Hapla   HDF5ReadCtx    h;
85969a06e7bSVaclav Hapla   PetscErrorCode ierr;
86069a06e7bSVaclav Hapla 
86169a06e7bSVaclav Hapla   PetscFunctionBegin;
86269a06e7bSVaclav Hapla   h = *ctx;
863b8ef406cSVaclav Hapla   PetscStackCallHDF5(H5Pclose,(h->plist));
86469a06e7bSVaclav Hapla   if (h->group != h->file) PetscStackCallHDF5(H5Gclose,(h->group));
86569a06e7bSVaclav Hapla   PetscStackCallHDF5(H5Sclose,(h->dataspace));
86669a06e7bSVaclav Hapla   PetscStackCallHDF5(H5Dclose,(h->dataset));
86769a06e7bSVaclav Hapla   ierr = PetscFree(*ctx);CHKERRQ(ierr);
86869a06e7bSVaclav Hapla   PetscFunctionReturn(0);
86969a06e7bSVaclav Hapla }
87069a06e7bSVaclav Hapla 
871eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadSizes_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout *map_)
87269a06e7bSVaclav Hapla {
87369a06e7bSVaclav Hapla   int            rdim, dim;
87469a06e7bSVaclav Hapla   hsize_t        dims[4];
87545e21b7fSVaclav Hapla   PetscInt       bsInd, lenInd, bs, N;
8768374c777SVaclav Hapla   PetscLayout    map;
8778374c777SVaclav Hapla   PetscErrorCode ierr;
87869a06e7bSVaclav Hapla 
87969a06e7bSVaclav Hapla   PetscFunctionBegin;
8808374c777SVaclav Hapla   if (!(*map_)) {
8818374c777SVaclav Hapla     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)viewer),map_);CHKERRQ(ierr);
8828374c777SVaclav Hapla   }
8838374c777SVaclav Hapla   map = *map_;
8849787f08bSVaclav Hapla   /* calculate expected number of dimensions */
885a5e1feadSVaclav Hapla   dim = 0;
8869568d583SVaclav Hapla   if (ctx->timestep >= 0) ++dim;
887a5e1feadSVaclav Hapla   ++dim; /* length in blocks */
8889568d583SVaclav Hapla   if (ctx->complexVal) ++dim;
8899787f08bSVaclav Hapla   /* get actual number of dimensions in dataset */
89069a06e7bSVaclav Hapla   PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(ctx->dataspace, dims, NULL));
8919787f08bSVaclav Hapla   /* calculate expected dimension indices */
8929787f08bSVaclav Hapla   lenInd = 0;
8939568d583SVaclav Hapla   if (ctx->timestep >= 0) ++lenInd;
8949787f08bSVaclav Hapla   bsInd = lenInd + 1;
8959568d583SVaclav Hapla   ctx->dim2 = PETSC_FALSE;
8969787f08bSVaclav Hapla   if (rdim == dim) {
89745e21b7fSVaclav Hapla     bs = 1; /* support vectors stored as 1D array */
8989787f08bSVaclav Hapla   } else if (rdim == dim+1) {
89945e21b7fSVaclav Hapla     bs = (PetscInt) dims[bsInd];
9009568d583SVaclav Hapla     if (bs == 1) ctx->dim2 = PETSC_TRUE; /* vector with blocksize of 1, still stored as 2D array */
901a5e1feadSVaclav Hapla   } else {
9029787f08bSVaclav Hapla     SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Dimension of array in file %d not %d as expected", rdim, dim);
903a5e1feadSVaclav Hapla   }
90445e21b7fSVaclav Hapla   N = (PetscInt) dims[lenInd]*bs;
9058374c777SVaclav Hapla 
9068374c777SVaclav Hapla   /* Set Vec sizes,blocksize,and type if not already set */
9078374c777SVaclav Hapla   if (map->bs < 0) {
90845e21b7fSVaclav Hapla     ierr = PetscLayoutSetBlockSize(map, bs);CHKERRQ(ierr);
90945e21b7fSVaclav 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);
9108374c777SVaclav Hapla   if (map->N < 0) {
91145e21b7fSVaclav Hapla     ierr = PetscLayoutSetSize(map, N);CHKERRQ(ierr);
91245e21b7fSVaclav 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);
91369a06e7bSVaclav Hapla   PetscFunctionReturn(0);
91469a06e7bSVaclav Hapla }
91569a06e7bSVaclav Hapla 
916eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadSelectHyperslab_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout map, hid_t *memspace)
91716f6402dSVaclav Hapla {
91816f6402dSVaclav Hapla   hsize_t        count[4], offset[4];
91916f6402dSVaclav Hapla   int            dim;
92016f6402dSVaclav Hapla   PetscInt       bs, n, low;
92116f6402dSVaclav Hapla   PetscErrorCode ierr;
92216f6402dSVaclav Hapla 
92316f6402dSVaclav Hapla   PetscFunctionBegin;
92416f6402dSVaclav Hapla   /* Compute local size and ownership range */
92516f6402dSVaclav Hapla   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
92616f6402dSVaclav Hapla   ierr = PetscLayoutGetBlockSize(map, &bs);CHKERRQ(ierr);
92716f6402dSVaclav Hapla   ierr = PetscLayoutGetLocalSize(map, &n);CHKERRQ(ierr);
92816f6402dSVaclav Hapla   ierr = PetscLayoutGetRange(map, &low, NULL);CHKERRQ(ierr);
92916f6402dSVaclav Hapla 
93016f6402dSVaclav Hapla   /* Each process defines a dataset and reads it from the hyperslab in the file */
93116f6402dSVaclav Hapla   dim  = 0;
93216f6402dSVaclav Hapla   if (ctx->timestep >= 0) {
93316f6402dSVaclav Hapla     count[dim]  = 1;
93416f6402dSVaclav Hapla     offset[dim] = ctx->timestep;
93516f6402dSVaclav Hapla     ++dim;
93616f6402dSVaclav Hapla   }
93716f6402dSVaclav Hapla   {
93816f6402dSVaclav Hapla     ierr = PetscHDF5IntCast(n/bs, &count[dim]);CHKERRQ(ierr);
93916f6402dSVaclav Hapla     ierr = PetscHDF5IntCast(low/bs, &offset[dim]);CHKERRQ(ierr);
94016f6402dSVaclav Hapla     ++dim;
94116f6402dSVaclav Hapla   }
94216f6402dSVaclav Hapla   if (bs > 1 || ctx->dim2) {
94316f6402dSVaclav Hapla     count[dim]  = bs;
94416f6402dSVaclav Hapla     offset[dim] = 0;
94516f6402dSVaclav Hapla     ++dim;
94616f6402dSVaclav Hapla   }
94716f6402dSVaclav Hapla   if (ctx->complexVal) {
94816f6402dSVaclav Hapla     count[dim]  = 2;
94916f6402dSVaclav Hapla     offset[dim] = 0;
95016f6402dSVaclav Hapla     ++dim;
95116f6402dSVaclav Hapla   }
95216f6402dSVaclav Hapla   PetscStackCallHDF5Return(*memspace,H5Screate_simple,(dim, count, NULL));
95316f6402dSVaclav Hapla   PetscStackCallHDF5(H5Sselect_hyperslab,(ctx->dataspace, H5S_SELECT_SET, offset, NULL, count, NULL));
95416f6402dSVaclav Hapla   PetscFunctionReturn(0);
95516f6402dSVaclav Hapla }
95616f6402dSVaclav Hapla 
957eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadArray_Private(PetscViewer viewer, HDF5ReadCtx h, hid_t datatype, hid_t memspace, void *arr)
958ef2d82ceSVaclav Hapla {
959ef2d82ceSVaclav Hapla   PetscFunctionBegin;
960ef2d82ceSVaclav Hapla   PetscStackCallHDF5(H5Dread,(h->dataset, datatype, memspace, h->dataspace, h->plist, arr));
961ef2d82ceSVaclav Hapla   PetscFunctionReturn(0);
962ef2d82ceSVaclav Hapla }
963ef2d82ceSVaclav Hapla 
964*dad982a8SVaclav Hapla PetscErrorCode PetscViewerHDF5Load(PetscViewer viewer, const char *name, PetscLayout map, hid_t datatype, void **newarr)
965a25c73c6SVaclav Hapla {
966a25c73c6SVaclav Hapla   HDF5ReadCtx     h;
967a25c73c6SVaclav Hapla   hid_t           memspace;
968a25c73c6SVaclav Hapla   size_t          unitsize;
969a25c73c6SVaclav Hapla   void            *arr;
970a25c73c6SVaclav Hapla   PetscErrorCode  ierr;
971a25c73c6SVaclav Hapla 
972a25c73c6SVaclav Hapla   PetscFunctionBegin;
973eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr);
974eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr);
975a25c73c6SVaclav Hapla   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
976eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadSelectHyperslab_Private(viewer, h, map, &memspace);CHKERRQ(ierr);
9774fc17bcdSVaclav Hapla 
9785a89fdf4SVaclav Hapla #if defined(PETSC_USE_COMPLEX)
9795a89fdf4SVaclav Hapla   if (!h->complexVal) {
980c76551afSVaclav Hapla     H5T_class_t clazz = H5Tget_class(datatype);
981c76551afSVaclav 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.");
9825a89fdf4SVaclav Hapla   }
9835a89fdf4SVaclav Hapla #else
9845a89fdf4SVaclav 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.");
9855a89fdf4SVaclav Hapla #endif
9864fc17bcdSVaclav Hapla   unitsize = H5Tget_size(datatype);
9874fc17bcdSVaclav Hapla   if (h->complexVal) unitsize *= 2;
9884fc17bcdSVaclav Hapla   ierr = PetscMalloc(map->n*unitsize, &arr);CHKERRQ(ierr);
9894fc17bcdSVaclav Hapla 
990eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadArray_Private(viewer, h, datatype, memspace, arr);CHKERRQ(ierr);
991a25c73c6SVaclav Hapla   PetscStackCallHDF5(H5Sclose,(memspace));
992eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr);
993a25c73c6SVaclav Hapla   *newarr = arr;
994a25c73c6SVaclav Hapla   PetscFunctionReturn(0);
995a25c73c6SVaclav Hapla }
996a25c73c6SVaclav Hapla 
997c1aaad9cSVaclav Hapla /*@C
998c1aaad9cSVaclav Hapla  PetscViewerHDF5ReadSizes - Read block size and global size of a vector (Vec or IS) stored in an HDF5 file.
999c1aaad9cSVaclav Hapla 
1000c1aaad9cSVaclav Hapla   Input Parameters:
1001c1aaad9cSVaclav Hapla + viewer - The HDF5 viewer
1002c1aaad9cSVaclav Hapla - name   - The vector name
1003c1aaad9cSVaclav Hapla 
1004c1aaad9cSVaclav Hapla   Output Parameter:
1005c1aaad9cSVaclav Hapla + bs     - block size
1006c1aaad9cSVaclav Hapla - N      - global size
1007c1aaad9cSVaclav Hapla 
1008c1aaad9cSVaclav Hapla   Note:
1009c1aaad9cSVaclav Hapla   A vector is stored as an HDF5 dataspace with 1-4 dimensions in this order:
1010c1aaad9cSVaclav Hapla   1) # timesteps (optional), 2) # blocks, 3) # elements per block (optional), 4) real and imaginary part (only for complex).
1011c1aaad9cSVaclav Hapla 
1012c1aaad9cSVaclav Hapla   A vectors can be stored as a 2D dataspace even if its blocksize is 1; see PetscViewerHDF5SetBaseDimension2().
1013c1aaad9cSVaclav Hapla 
1014c1aaad9cSVaclav Hapla   Level: advanced
1015c1aaad9cSVaclav Hapla 
1016c1aaad9cSVaclav Hapla .seealso: PetscViewerHDF5Open(), VecLoad(), ISLoad(), VecGetSize(), ISGetSize(), PetscViewerHDF5SetBaseDimension2()
1017c1aaad9cSVaclav Hapla @*/
101869a06e7bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadSizes(PetscViewer viewer, const char name[], PetscInt *bs, PetscInt *N)
101969a06e7bSVaclav Hapla {
102069a06e7bSVaclav Hapla   HDF5ReadCtx    h;
10218374c777SVaclav Hapla   PetscLayout    map=NULL;
102269a06e7bSVaclav Hapla   PetscErrorCode ierr;
102369a06e7bSVaclav Hapla 
102469a06e7bSVaclav Hapla   PetscFunctionBegin;
1025c1aaad9cSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1026eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr);
1027eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr);
1028eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr);
10298374c777SVaclav Hapla   if (bs) *bs = map->bs;
10308374c777SVaclav Hapla   if (N) *N = map->N;
10318374c777SVaclav Hapla   ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr);
1032a5e1feadSVaclav Hapla   PetscFunctionReturn(0);
1033a5e1feadSVaclav Hapla }
1034a5e1feadSVaclav Hapla 
1035a75e6a4aSMatthew G. Knepley /*
1036a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1037a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
1038a75e6a4aSMatthew G. Knepley */
1039d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1040a75e6a4aSMatthew G. Knepley 
1041a75e6a4aSMatthew G. Knepley /*@C
1042a75e6a4aSMatthew G. Knepley   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
1043a75e6a4aSMatthew G. Knepley 
1044a75e6a4aSMatthew G. Knepley   Collective on MPI_Comm
1045a75e6a4aSMatthew G. Knepley 
1046a75e6a4aSMatthew G. Knepley   Input Parameter:
1047a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer
1048a75e6a4aSMatthew G. Knepley 
1049a75e6a4aSMatthew G. Knepley   Level: intermediate
1050a75e6a4aSMatthew G. Knepley 
1051a75e6a4aSMatthew G. Knepley   Options Database Keys:
1052a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name>
1053a75e6a4aSMatthew G. Knepley 
1054a75e6a4aSMatthew G. Knepley   Environmental variables:
1055a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME
1056a75e6a4aSMatthew G. Knepley 
1057a75e6a4aSMatthew G. Knepley   Notes:
1058a75e6a4aSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
1059a75e6a4aSMatthew G. Knepley   an error code.  The HDF5 PetscViewer is usually used in the form
1060a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1061a75e6a4aSMatthew G. Knepley 
1062a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1063a75e6a4aSMatthew G. Knepley @*/
1064a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1065a75e6a4aSMatthew G. Knepley {
1066a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
1067a75e6a4aSMatthew G. Knepley   PetscBool      flg;
1068a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
1069a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
1070a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
1071a75e6a4aSMatthew G. Knepley 
1072a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
1073a75e6a4aSMatthew 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);}
1074a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
107512801b39SBarry Smith     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
1076a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1077a75e6a4aSMatthew G. Knepley   }
107847435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1079a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1080a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
1081a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
1082a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1083a75e6a4aSMatthew G. Knepley     if (!flg) {
1084a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname,"output.h5");
1085a75e6a4aSMatthew G. Knepley       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1086a75e6a4aSMatthew G. Knepley     }
1087a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
1088a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1089a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
1090a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
109147435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
1092a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1093a75e6a4aSMatthew G. Knepley   }
1094a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
1095a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1096a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
1097a75e6a4aSMatthew G. Knepley }
1098