xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 45e21b7fbc366ef544fddda1e2bee75e89090bb7)
1af0996ceSBarry Smith #include <petsc/private/viewerimpl.h>    /*I   "petscsys.h"   I*/
2d70abbfaSBarry Smith #include <petscviewerhdf5.h>    /*I   "petscviewerhdf5.h"   I*/
35c6c1daeSBarry Smith 
45c6c1daeSBarry Smith typedef struct GroupList {
55c6c1daeSBarry Smith   const char       *name;
65c6c1daeSBarry Smith   struct GroupList *next;
75c6c1daeSBarry Smith } GroupList;
85c6c1daeSBarry Smith 
95c6c1daeSBarry Smith typedef struct {
105c6c1daeSBarry Smith   char          *filename;
115c6c1daeSBarry Smith   PetscFileMode btype;
125c6c1daeSBarry Smith   hid_t         file_id;
135c6c1daeSBarry Smith   PetscInt      timestep;
145c6c1daeSBarry Smith   GroupList     *groups;
1582ea9b62SBarry Smith   PetscBool     basedimension2;  /* save vectors and DMDA vectors with a dimension of at least 2 even if the bs/dof is 1 */
169a0502c6SHåkon Strandenes   PetscBool     spoutput;  /* write data in single precision even if PETSc is compiled with double precision PetscReal */
175c6c1daeSBarry Smith } PetscViewer_HDF5;
185c6c1daeSBarry Smith 
194416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v)
2082ea9b62SBarry Smith {
2182ea9b62SBarry Smith   PetscErrorCode   ierr;
2282ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;
2382ea9b62SBarry Smith 
2482ea9b62SBarry Smith   PetscFunctionBegin;
2582ea9b62SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr);
2682ea9b62SBarry Smith   ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr);
279a0502c6SHåkon Strandenes   ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr);
2882ea9b62SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
2982ea9b62SBarry Smith   PetscFunctionReturn(0);
3082ea9b62SBarry Smith }
3182ea9b62SBarry Smith 
325c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
335c6c1daeSBarry Smith {
345c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
355c6c1daeSBarry Smith   PetscErrorCode   ierr;
365c6c1daeSBarry Smith 
375c6c1daeSBarry Smith   PetscFunctionBegin;
385c6c1daeSBarry Smith   ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);
39729ab6d9SBarry Smith   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
405c6c1daeSBarry Smith   PetscFunctionReturn(0);
415c6c1daeSBarry Smith }
425c6c1daeSBarry Smith 
435c6c1daeSBarry Smith PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
445c6c1daeSBarry Smith {
455c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
465c6c1daeSBarry Smith   PetscErrorCode   ierr;
475c6c1daeSBarry Smith 
485c6c1daeSBarry Smith   PetscFunctionBegin;
495c6c1daeSBarry Smith   ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr);
505c6c1daeSBarry Smith   while (hdf5->groups) {
515c6c1daeSBarry Smith     GroupList *tmp = hdf5->groups->next;
525c6c1daeSBarry Smith 
535c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups->name);CHKERRQ(ierr);
545c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups);CHKERRQ(ierr);
555c6c1daeSBarry Smith     hdf5->groups = tmp;
565c6c1daeSBarry Smith   }
575c6c1daeSBarry Smith   ierr = PetscFree(hdf5);CHKERRQ(ierr);
580b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
59d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr);
600b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
615c6c1daeSBarry Smith   PetscFunctionReturn(0);
625c6c1daeSBarry Smith }
635c6c1daeSBarry Smith 
645c6c1daeSBarry Smith PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
655c6c1daeSBarry Smith {
665c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
675c6c1daeSBarry Smith 
685c6c1daeSBarry Smith   PetscFunctionBegin;
695c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
705c6c1daeSBarry Smith   hdf5->btype = type;
715c6c1daeSBarry Smith   PetscFunctionReturn(0);
725c6c1daeSBarry Smith }
735c6c1daeSBarry Smith 
7482ea9b62SBarry Smith PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
7582ea9b62SBarry Smith {
7682ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7782ea9b62SBarry Smith 
7882ea9b62SBarry Smith   PetscFunctionBegin;
7982ea9b62SBarry Smith   hdf5->basedimension2 = flg;
8082ea9b62SBarry Smith   PetscFunctionReturn(0);
8182ea9b62SBarry Smith }
8282ea9b62SBarry Smith 
83df863907SAlex Fikl /*@
8482ea9b62SBarry Smith      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
8582ea9b62SBarry Smith        dimension of 2.
8682ea9b62SBarry Smith 
8782ea9b62SBarry Smith     Logically Collective on PetscViewer
8882ea9b62SBarry Smith 
8982ea9b62SBarry Smith   Input Parameters:
9082ea9b62SBarry Smith +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
9182ea9b62SBarry Smith -  flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1
9282ea9b62SBarry Smith 
9382ea9b62SBarry Smith   Options Database:
9482ea9b62SBarry Smith .  -viewer_hdf5_base_dimension2 - turns on (true) or off (false) using a dimension of 2 in the HDF5 file even if the bs/dof of the vector is 1
9582ea9b62SBarry Smith 
9682ea9b62SBarry Smith 
9795452b02SPatrick Sanan   Notes:
9895452b02SPatrick 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
9982ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
10082ea9b62SBarry Smith 
10182ea9b62SBarry Smith   Level: intermediate
10282ea9b62SBarry Smith 
10382ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
10482ea9b62SBarry Smith 
10582ea9b62SBarry Smith @*/
10682ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
10782ea9b62SBarry Smith {
10882ea9b62SBarry Smith   PetscErrorCode ierr;
10982ea9b62SBarry Smith 
11082ea9b62SBarry Smith   PetscFunctionBegin;
11182ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
11282ea9b62SBarry Smith   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
11382ea9b62SBarry Smith   PetscFunctionReturn(0);
11482ea9b62SBarry Smith }
11582ea9b62SBarry Smith 
116df863907SAlex Fikl /*@
11782ea9b62SBarry Smith      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
11882ea9b62SBarry Smith        dimension of 2.
11982ea9b62SBarry Smith 
12082ea9b62SBarry Smith     Logically Collective on PetscViewer
12182ea9b62SBarry Smith 
12282ea9b62SBarry Smith   Input Parameter:
12382ea9b62SBarry Smith .  viewer - the PetscViewer, must be of type HDF5
12482ea9b62SBarry Smith 
12582ea9b62SBarry Smith   Output Parameter:
12682ea9b62SBarry 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
12782ea9b62SBarry Smith 
12895452b02SPatrick Sanan   Notes:
12995452b02SPatrick 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
13082ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
13182ea9b62SBarry Smith 
13282ea9b62SBarry Smith   Level: intermediate
13382ea9b62SBarry Smith 
13482ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
13582ea9b62SBarry Smith 
13682ea9b62SBarry Smith @*/
13782ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
13882ea9b62SBarry Smith {
13982ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
14082ea9b62SBarry Smith 
14182ea9b62SBarry Smith   PetscFunctionBegin;
14282ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
14382ea9b62SBarry Smith   *flg = hdf5->basedimension2;
14482ea9b62SBarry Smith   PetscFunctionReturn(0);
14582ea9b62SBarry Smith }
14682ea9b62SBarry Smith 
1479a0502c6SHåkon Strandenes PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
1489a0502c6SHåkon Strandenes {
1499a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1509a0502c6SHåkon Strandenes 
1519a0502c6SHåkon Strandenes   PetscFunctionBegin;
1529a0502c6SHåkon Strandenes   hdf5->spoutput = flg;
1539a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
1549a0502c6SHåkon Strandenes }
1559a0502c6SHåkon Strandenes 
156df863907SAlex Fikl /*@
1579a0502c6SHåkon Strandenes      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
1589a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
1599a0502c6SHåkon Strandenes 
1609a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
1619a0502c6SHåkon Strandenes 
1629a0502c6SHåkon Strandenes   Input Parameters:
1639a0502c6SHåkon Strandenes +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
1649a0502c6SHåkon Strandenes -  flg - if PETSC_TRUE the data will be written to disk with single precision
1659a0502c6SHåkon Strandenes 
1669a0502c6SHåkon Strandenes   Options Database:
1679a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
1689a0502c6SHåkon Strandenes 
1699a0502c6SHåkon Strandenes 
17095452b02SPatrick Sanan   Notes:
17195452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
1729a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
1739a0502c6SHåkon Strandenes 
1749a0502c6SHåkon Strandenes   Level: intermediate
1759a0502c6SHåkon Strandenes 
1769a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
1779a0502c6SHåkon Strandenes           PetscReal
1789a0502c6SHåkon Strandenes 
1799a0502c6SHåkon Strandenes @*/
1809a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
1819a0502c6SHåkon Strandenes {
1829a0502c6SHåkon Strandenes   PetscErrorCode ierr;
1839a0502c6SHåkon Strandenes 
1849a0502c6SHåkon Strandenes   PetscFunctionBegin;
1859a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1869a0502c6SHåkon Strandenes   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
1879a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
1889a0502c6SHåkon Strandenes }
1899a0502c6SHåkon Strandenes 
190df863907SAlex Fikl /*@
1919a0502c6SHåkon Strandenes      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
1929a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
1939a0502c6SHåkon Strandenes 
1949a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
1959a0502c6SHåkon Strandenes 
1969a0502c6SHåkon Strandenes   Input Parameter:
1979a0502c6SHåkon Strandenes .  viewer - the PetscViewer, must be of type HDF5
1989a0502c6SHåkon Strandenes 
1999a0502c6SHåkon Strandenes   Output Parameter:
2009a0502c6SHåkon Strandenes .  flg - if PETSC_TRUE the data will be written to disk with single precision
2019a0502c6SHåkon Strandenes 
20295452b02SPatrick Sanan   Notes:
20395452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2049a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2059a0502c6SHåkon Strandenes 
2069a0502c6SHåkon Strandenes   Level: intermediate
2079a0502c6SHåkon Strandenes 
2089a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2099a0502c6SHåkon Strandenes           PetscReal
2109a0502c6SHåkon Strandenes 
2119a0502c6SHåkon Strandenes @*/
2129a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
2139a0502c6SHåkon Strandenes {
2149a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2159a0502c6SHåkon Strandenes 
2169a0502c6SHåkon Strandenes   PetscFunctionBegin;
2179a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2189a0502c6SHåkon Strandenes   *flg = hdf5->spoutput;
2199a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2209a0502c6SHåkon Strandenes }
2219a0502c6SHåkon Strandenes 
2225c6c1daeSBarry Smith PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
2235c6c1daeSBarry Smith {
2245c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2255c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
2265c6c1daeSBarry Smith   MPI_Info          info = MPI_INFO_NULL;
2275c6c1daeSBarry Smith #endif
2285c6c1daeSBarry Smith   hid_t             plist_id;
2295c6c1daeSBarry Smith   PetscErrorCode    ierr;
2305c6c1daeSBarry Smith 
2315c6c1daeSBarry Smith   PetscFunctionBegin;
232f683cc58SToby Isaac   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
233f683cc58SToby Isaac   if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);}
2345c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
2355c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
236729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
2375c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
238729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info));
2395c6c1daeSBarry Smith #endif
2405c6c1daeSBarry Smith   /* Create or open the file collectively */
2415c6c1daeSBarry Smith   switch (hdf5->btype) {
2425c6c1daeSBarry Smith   case FILE_MODE_READ:
243729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
2445c6c1daeSBarry Smith     break;
2455c6c1daeSBarry Smith   case FILE_MODE_APPEND:
246729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
2475c6c1daeSBarry Smith     break;
2485c6c1daeSBarry Smith   case FILE_MODE_WRITE:
249729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
2505c6c1daeSBarry Smith     break;
2515c6c1daeSBarry Smith   default:
2525c6c1daeSBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
2535c6c1daeSBarry Smith   }
2545c6c1daeSBarry Smith   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
255729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
2565c6c1daeSBarry Smith   PetscFunctionReturn(0);
2575c6c1daeSBarry Smith }
2585c6c1daeSBarry Smith 
259d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
260d1232d7fSToby Isaac {
261d1232d7fSToby Isaac   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
262d1232d7fSToby Isaac 
263d1232d7fSToby Isaac   PetscFunctionBegin;
264d1232d7fSToby Isaac   *name = vhdf5->filename;
265d1232d7fSToby Isaac   PetscFunctionReturn(0);
266d1232d7fSToby Isaac }
267d1232d7fSToby Isaac 
2688556b5ebSBarry Smith /*MC
2698556b5ebSBarry Smith    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
2708556b5ebSBarry Smith 
2718556b5ebSBarry Smith 
2728556b5ebSBarry Smith .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
2738556b5ebSBarry Smith            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
2748556b5ebSBarry Smith            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
2758556b5ebSBarry Smith            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
2768556b5ebSBarry Smith 
2771b266c99SBarry Smith   Level: beginner
2788556b5ebSBarry Smith M*/
279d1232d7fSToby Isaac 
2808cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
2815c6c1daeSBarry Smith {
2825c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
2835c6c1daeSBarry Smith   PetscErrorCode   ierr;
2845c6c1daeSBarry Smith 
2855c6c1daeSBarry Smith   PetscFunctionBegin;
286b00a9115SJed Brown   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
2875c6c1daeSBarry Smith 
2885c6c1daeSBarry Smith   v->data                = (void*) hdf5;
2895c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
29082ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
2915c6c1daeSBarry Smith   v->ops->flush          = 0;
2925c6c1daeSBarry Smith   hdf5->btype            = (PetscFileMode) -1;
2935c6c1daeSBarry Smith   hdf5->filename         = 0;
2945c6c1daeSBarry Smith   hdf5->timestep         = -1;
2950298fd71SBarry Smith   hdf5->groups           = NULL;
2965c6c1daeSBarry Smith 
2970b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
298d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr);
2990b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
30082ea9b62SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
3019a0502c6SHåkon Strandenes   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr);
3025c6c1daeSBarry Smith   PetscFunctionReturn(0);
3035c6c1daeSBarry Smith }
3045c6c1daeSBarry Smith 
3055c6c1daeSBarry Smith /*@C
3065c6c1daeSBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
3075c6c1daeSBarry Smith 
3085c6c1daeSBarry Smith    Collective on MPI_Comm
3095c6c1daeSBarry Smith 
3105c6c1daeSBarry Smith    Input Parameters:
3115c6c1daeSBarry Smith +  comm - MPI communicator
3125c6c1daeSBarry Smith .  name - name of file
3135c6c1daeSBarry Smith -  type - type of file
3145c6c1daeSBarry Smith $    FILE_MODE_WRITE - create new file for binary output
3155c6c1daeSBarry Smith $    FILE_MODE_READ - open existing file for binary input
3165c6c1daeSBarry Smith $    FILE_MODE_APPEND - open existing file for binary output
3175c6c1daeSBarry Smith 
3185c6c1daeSBarry Smith    Output Parameter:
3195c6c1daeSBarry Smith .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
3205c6c1daeSBarry Smith 
32182ea9b62SBarry Smith   Options Database:
32282ea9b62SBarry 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
3239a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
32482ea9b62SBarry Smith 
3255c6c1daeSBarry Smith    Level: beginner
3265c6c1daeSBarry Smith 
3275c6c1daeSBarry Smith    Note:
3285c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
3295c6c1daeSBarry Smith 
3305c6c1daeSBarry Smith    Concepts: HDF5 files
3315c6c1daeSBarry Smith    Concepts: PetscViewerHDF5^creating
3325c6c1daeSBarry Smith 
3336a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
3349a0502c6SHåkon Strandenes           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
335a56f64adSBarry Smith           MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
3365c6c1daeSBarry Smith @*/
3375c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
3385c6c1daeSBarry Smith {
3395c6c1daeSBarry Smith   PetscErrorCode ierr;
3405c6c1daeSBarry Smith 
3415c6c1daeSBarry Smith   PetscFunctionBegin;
3425c6c1daeSBarry Smith   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
3435c6c1daeSBarry Smith   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
3445c6c1daeSBarry Smith   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
3455c6c1daeSBarry Smith   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
3465c6c1daeSBarry Smith   PetscFunctionReturn(0);
3475c6c1daeSBarry Smith }
3485c6c1daeSBarry Smith 
3495c6c1daeSBarry Smith /*@C
3505c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
3515c6c1daeSBarry Smith 
3525c6c1daeSBarry Smith   Not collective
3535c6c1daeSBarry Smith 
3545c6c1daeSBarry Smith   Input Parameter:
3555c6c1daeSBarry Smith . viewer - the PetscViewer
3565c6c1daeSBarry Smith 
3575c6c1daeSBarry Smith   Output Parameter:
3585c6c1daeSBarry Smith . file_id - The file id
3595c6c1daeSBarry Smith 
3605c6c1daeSBarry Smith   Level: intermediate
3615c6c1daeSBarry Smith 
3625c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open()
3635c6c1daeSBarry Smith @*/
3645c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
3655c6c1daeSBarry Smith {
3665c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3675c6c1daeSBarry Smith 
3685c6c1daeSBarry Smith   PetscFunctionBegin;
3695c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
3705c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
3715c6c1daeSBarry Smith   PetscFunctionReturn(0);
3725c6c1daeSBarry Smith }
3735c6c1daeSBarry Smith 
3745c6c1daeSBarry Smith /*@C
3755c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
3765c6c1daeSBarry Smith 
3775c6c1daeSBarry Smith   Not collective
3785c6c1daeSBarry Smith 
3795c6c1daeSBarry Smith   Input Parameters:
3805c6c1daeSBarry Smith + viewer - the PetscViewer
3815c6c1daeSBarry Smith - name - The group name
3825c6c1daeSBarry Smith 
3835c6c1daeSBarry Smith   Level: intermediate
3845c6c1daeSBarry Smith 
385874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
3865c6c1daeSBarry Smith @*/
3875c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
3885c6c1daeSBarry Smith {
3895c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3905c6c1daeSBarry Smith   GroupList        *groupNode;
3915c6c1daeSBarry Smith   PetscErrorCode   ierr;
3925c6c1daeSBarry Smith 
3935c6c1daeSBarry Smith   PetscFunctionBegin;
3945c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
3955c6c1daeSBarry Smith   PetscValidCharPointer(name,2);
39695dccacaSBarry Smith   ierr = PetscNew(&groupNode);CHKERRQ(ierr);
3975c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
398a297a907SKarl Rupp 
3995c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
4005c6c1daeSBarry Smith   hdf5->groups    = groupNode;
4015c6c1daeSBarry Smith   PetscFunctionReturn(0);
4025c6c1daeSBarry Smith }
4035c6c1daeSBarry Smith 
4043ef9c667SSatish Balay /*@
4055c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
4065c6c1daeSBarry Smith 
4075c6c1daeSBarry Smith   Not collective
4085c6c1daeSBarry Smith 
4095c6c1daeSBarry Smith   Input Parameter:
4105c6c1daeSBarry Smith . viewer - the PetscViewer
4115c6c1daeSBarry Smith 
4125c6c1daeSBarry Smith   Level: intermediate
4135c6c1daeSBarry Smith 
414874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
4155c6c1daeSBarry Smith @*/
4165c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
4175c6c1daeSBarry Smith {
4185c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4195c6c1daeSBarry Smith   GroupList        *groupNode;
4205c6c1daeSBarry Smith   PetscErrorCode   ierr;
4215c6c1daeSBarry Smith 
4225c6c1daeSBarry Smith   PetscFunctionBegin;
4235c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
42482f516ccSBarry Smith   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
4255c6c1daeSBarry Smith   groupNode    = hdf5->groups;
4265c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
4275c6c1daeSBarry Smith   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
4285c6c1daeSBarry Smith   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
4295c6c1daeSBarry Smith   PetscFunctionReturn(0);
4305c6c1daeSBarry Smith }
4315c6c1daeSBarry Smith 
4325c6c1daeSBarry Smith /*@C
433874270d9SVaclav Hapla   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup().
434874270d9SVaclav Hapla   If none has been assigned, returns NULL.
4355c6c1daeSBarry Smith 
4365c6c1daeSBarry Smith   Not collective
4375c6c1daeSBarry Smith 
4385c6c1daeSBarry Smith   Input Parameter:
4395c6c1daeSBarry Smith . viewer - the PetscViewer
4405c6c1daeSBarry Smith 
4415c6c1daeSBarry Smith   Output Parameter:
4425c6c1daeSBarry Smith . name - The group name
4435c6c1daeSBarry Smith 
4445c6c1daeSBarry Smith   Level: intermediate
4455c6c1daeSBarry Smith 
446874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
4475c6c1daeSBarry Smith @*/
4485c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
4495c6c1daeSBarry Smith {
4505c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
4515c6c1daeSBarry Smith 
4525c6c1daeSBarry Smith   PetscFunctionBegin;
4535c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
454c959eef4SJed Brown   PetscValidPointer(name,2);
455a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
4560298fd71SBarry Smith   else *name = NULL;
4575c6c1daeSBarry Smith   PetscFunctionReturn(0);
4585c6c1daeSBarry Smith }
4595c6c1daeSBarry Smith 
460874270d9SVaclav Hapla /*@C
461874270d9SVaclav Hapla   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(),
462874270d9SVaclav Hapla   and return this group's ID and file ID.
463874270d9SVaclav Hapla   If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID.
464874270d9SVaclav Hapla 
465874270d9SVaclav Hapla   Not collective
466874270d9SVaclav Hapla 
467874270d9SVaclav Hapla   Input Parameter:
468874270d9SVaclav Hapla . viewer - the PetscViewer
469874270d9SVaclav Hapla 
470874270d9SVaclav Hapla   Output Parameter:
471874270d9SVaclav Hapla + fileId - The HDF5 file ID
472874270d9SVaclav Hapla - groupId - The HDF5 group ID
473874270d9SVaclav Hapla 
474874270d9SVaclav Hapla   Level: intermediate
475874270d9SVaclav Hapla 
476874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
477874270d9SVaclav Hapla @*/
47854dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
47954dbf706SVaclav Hapla {
48054dbf706SVaclav Hapla   hid_t          file_id, group;
48154dbf706SVaclav Hapla   htri_t         found;
48254dbf706SVaclav Hapla   const char     *groupName = NULL;
48354dbf706SVaclav Hapla   PetscErrorCode ierr;
48454dbf706SVaclav Hapla 
48554dbf706SVaclav Hapla   PetscFunctionBegin;
48654dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr);
48754dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr);
48854dbf706SVaclav Hapla   /* Open group */
48954dbf706SVaclav Hapla   if (groupName) {
49054dbf706SVaclav Hapla     PetscBool root;
49154dbf706SVaclav Hapla 
49254dbf706SVaclav Hapla     ierr = PetscStrcmp(groupName, "/", &root);CHKERRQ(ierr);
49354dbf706SVaclav Hapla     PetscStackCall("H5Lexists",found = H5Lexists(file_id, groupName, H5P_DEFAULT));
49454dbf706SVaclav Hapla     if (!root && (found <= 0)) {
49554dbf706SVaclav Hapla #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
49654dbf706SVaclav Hapla       PetscStackCallHDF5Return(group,H5Gcreate2,(file_id, groupName, 0, H5P_DEFAULT, H5P_DEFAULT));
49754dbf706SVaclav Hapla #else /* deprecated HDF5 1.6 API */
49854dbf706SVaclav Hapla       PetscStackCallHDF5Return(group,H5Gcreate,(file_id, groupName, 0));
49954dbf706SVaclav Hapla #endif
50054dbf706SVaclav Hapla       PetscStackCallHDF5(H5Gclose,(group));
50154dbf706SVaclav Hapla     }
50254dbf706SVaclav Hapla #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
50354dbf706SVaclav Hapla     PetscStackCallHDF5Return(group,H5Gopen2,(file_id, groupName, H5P_DEFAULT));
50454dbf706SVaclav Hapla #else
50554dbf706SVaclav Hapla     PetscStackCallHDF5Return(group,H5Gopen,(file_id, groupName));
50654dbf706SVaclav Hapla #endif
50754dbf706SVaclav Hapla   } else group = file_id;
50854dbf706SVaclav Hapla 
50954dbf706SVaclav Hapla   *fileId  = file_id;
51054dbf706SVaclav Hapla   *groupId = group;
51154dbf706SVaclav Hapla   PetscFunctionReturn(0);
51254dbf706SVaclav Hapla }
51354dbf706SVaclav Hapla 
5143ef9c667SSatish Balay /*@
5155c6c1daeSBarry Smith   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
5165c6c1daeSBarry Smith 
5175c6c1daeSBarry Smith   Not collective
5185c6c1daeSBarry Smith 
5195c6c1daeSBarry Smith   Input Parameter:
5205c6c1daeSBarry Smith . viewer - the PetscViewer
5215c6c1daeSBarry Smith 
5225c6c1daeSBarry Smith   Level: intermediate
5235c6c1daeSBarry Smith 
5245c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
5255c6c1daeSBarry Smith @*/
5265c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
5275c6c1daeSBarry Smith {
5285c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5295c6c1daeSBarry Smith 
5305c6c1daeSBarry Smith   PetscFunctionBegin;
5315c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5325c6c1daeSBarry Smith   ++hdf5->timestep;
5335c6c1daeSBarry Smith   PetscFunctionReturn(0);
5345c6c1daeSBarry Smith }
5355c6c1daeSBarry Smith 
5363ef9c667SSatish Balay /*@
5375c6c1daeSBarry Smith   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
5385c6c1daeSBarry Smith   of -1 disables blocking with timesteps.
5395c6c1daeSBarry Smith 
5405c6c1daeSBarry Smith   Not collective
5415c6c1daeSBarry Smith 
5425c6c1daeSBarry Smith   Input Parameters:
5435c6c1daeSBarry Smith + viewer - the PetscViewer
5445c6c1daeSBarry Smith - timestep - The timestep number
5455c6c1daeSBarry Smith 
5465c6c1daeSBarry Smith   Level: intermediate
5475c6c1daeSBarry Smith 
5485c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
5495c6c1daeSBarry Smith @*/
5505c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
5515c6c1daeSBarry Smith {
5525c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5535c6c1daeSBarry Smith 
5545c6c1daeSBarry Smith   PetscFunctionBegin;
5555c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5565c6c1daeSBarry Smith   hdf5->timestep = timestep;
5575c6c1daeSBarry Smith   PetscFunctionReturn(0);
5585c6c1daeSBarry Smith }
5595c6c1daeSBarry Smith 
5603ef9c667SSatish Balay /*@
5615c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
5625c6c1daeSBarry Smith 
5635c6c1daeSBarry Smith   Not collective
5645c6c1daeSBarry Smith 
5655c6c1daeSBarry Smith   Input Parameter:
5665c6c1daeSBarry Smith . viewer - the PetscViewer
5675c6c1daeSBarry Smith 
5685c6c1daeSBarry Smith   Output Parameter:
5695c6c1daeSBarry Smith . timestep - The timestep number
5705c6c1daeSBarry Smith 
5715c6c1daeSBarry Smith   Level: intermediate
5725c6c1daeSBarry Smith 
5735c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
5745c6c1daeSBarry Smith @*/
5755c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
5765c6c1daeSBarry Smith {
5775c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5785c6c1daeSBarry Smith 
5795c6c1daeSBarry Smith   PetscFunctionBegin;
5805c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5815c6c1daeSBarry Smith   PetscValidPointer(timestep,2);
5825c6c1daeSBarry Smith   *timestep = hdf5->timestep;
5835c6c1daeSBarry Smith   PetscFunctionReturn(0);
5845c6c1daeSBarry Smith }
5855c6c1daeSBarry Smith 
58636bce6e8SMatthew G. Knepley /*@C
58736bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
58836bce6e8SMatthew G. Knepley 
58936bce6e8SMatthew G. Knepley   Not collective
59036bce6e8SMatthew G. Knepley 
59136bce6e8SMatthew G. Knepley   Input Parameter:
59236bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
59336bce6e8SMatthew G. Knepley 
59436bce6e8SMatthew G. Knepley   Output Parameter:
59536bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
59636bce6e8SMatthew G. Knepley 
59736bce6e8SMatthew G. Knepley   Level: advanced
59836bce6e8SMatthew G. Knepley 
59936bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
60036bce6e8SMatthew G. Knepley @*/
60136bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
60236bce6e8SMatthew G. Knepley {
60336bce6e8SMatthew G. Knepley   PetscFunctionBegin;
60436bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
60536bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
60636bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_LLONG;
60736bce6e8SMatthew G. Knepley #else
60836bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_INT;
60936bce6e8SMatthew G. Knepley #endif
61036bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
61136bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
61236bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
61336bce6e8SMatthew G. Knepley   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_DOUBLE;
61436bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_DOUBLE;
61536bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
61636bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
61736bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
6187e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
61936bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
62036bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
62136bce6e8SMatthew G. Knepley }
62236bce6e8SMatthew G. Knepley 
62336bce6e8SMatthew G. Knepley /*@C
62436bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
62536bce6e8SMatthew G. Knepley 
62636bce6e8SMatthew G. Knepley   Not collective
62736bce6e8SMatthew G. Knepley 
62836bce6e8SMatthew G. Knepley   Input Parameter:
62936bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
63036bce6e8SMatthew G. Knepley 
63136bce6e8SMatthew G. Knepley   Output Parameter:
63236bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
63336bce6e8SMatthew G. Knepley 
63436bce6e8SMatthew G. Knepley   Level: advanced
63536bce6e8SMatthew G. Knepley 
63636bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
63736bce6e8SMatthew G. Knepley @*/
63836bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
63936bce6e8SMatthew G. Knepley {
64036bce6e8SMatthew G. Knepley   PetscFunctionBegin;
64136bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
64236bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
64336bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
64436bce6e8SMatthew G. Knepley #else
64536bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
64636bce6e8SMatthew G. Knepley #endif
64736bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
64836bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
64936bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
65036bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
65136bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
65236bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
6537e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
65436bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
65536bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
65636bce6e8SMatthew G. Knepley }
65736bce6e8SMatthew G. Knepley 
658df863907SAlex Fikl /*@C
65936bce6e8SMatthew G. Knepley  PetscViewerHDF5WriteAttribute - Write a scalar attribute
66036bce6e8SMatthew G. Knepley 
66136bce6e8SMatthew G. Knepley   Input Parameters:
66236bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer
66336bce6e8SMatthew G. Knepley . parent - The parent name
66436bce6e8SMatthew G. Knepley . name   - The attribute name
66536bce6e8SMatthew G. Knepley . datatype - The attribute type
66636bce6e8SMatthew G. Knepley - value    - The attribute value
66736bce6e8SMatthew G. Knepley 
66836bce6e8SMatthew G. Knepley   Level: advanced
66936bce6e8SMatthew G. Knepley 
670e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute()
67136bce6e8SMatthew G. Knepley @*/
67236bce6e8SMatthew G. Knepley PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
67336bce6e8SMatthew G. Knepley {
67460568592SMatthew G. Knepley   hid_t          h5, dataspace, obj, attribute, dtype;
67536bce6e8SMatthew G. Knepley   PetscErrorCode ierr;
67636bce6e8SMatthew G. Knepley 
67736bce6e8SMatthew G. Knepley   PetscFunctionBegin;
6785cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
67936bce6e8SMatthew G. Knepley   PetscValidPointer(parent, 2);
68036bce6e8SMatthew G. Knepley   PetscValidPointer(name, 3);
68136bce6e8SMatthew G. Knepley   PetscValidPointer(value, 4);
68236bce6e8SMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
6837e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
6847e97c476SMatthew G. Knepley     size_t len;
6857e97c476SMatthew G. Knepley     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
686729ab6d9SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
6877e97c476SMatthew G. Knepley   }
68836bce6e8SMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
689729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
69060568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
69136bce6e8SMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
69260568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
69336bce6e8SMatthew G. Knepley #else
69460568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Acreate,(obj, name, dtype, dataspace, H5P_DEFAULT));
69536bce6e8SMatthew G. Knepley #endif
696729ab6d9SBarry Smith   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
697729ab6d9SBarry Smith   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
698729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
69960568592SMatthew G. Knepley   PetscStackCallHDF5(H5Oclose,(obj));
700729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(dataspace));
70136bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
70236bce6e8SMatthew G. Knepley }
70336bce6e8SMatthew G. Knepley 
704df863907SAlex Fikl /*@C
705ad0c61feSMatthew G. Knepley  PetscViewerHDF5ReadAttribute - Read a scalar attribute
706ad0c61feSMatthew G. Knepley 
707ad0c61feSMatthew G. Knepley   Input Parameters:
708ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer
709ad0c61feSMatthew G. Knepley . parent - The parent name
710ad0c61feSMatthew G. Knepley . name   - The attribute name
711ad0c61feSMatthew G. Knepley - datatype - The attribute type
712ad0c61feSMatthew G. Knepley 
713ad0c61feSMatthew G. Knepley   Output Parameter:
714ad0c61feSMatthew G. Knepley . value    - The attribute value
715ad0c61feSMatthew G. Knepley 
716ad0c61feSMatthew G. Knepley   Level: advanced
717ad0c61feSMatthew G. Knepley 
718e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute()
719ad0c61feSMatthew G. Knepley @*/
720ad0c61feSMatthew G. Knepley PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, void *value)
721ad0c61feSMatthew G. Knepley {
722f0b58479SMatthew G. Knepley   hid_t          h5, obj, attribute, atype, dtype;
723ad0c61feSMatthew G. Knepley   PetscErrorCode ierr;
724ad0c61feSMatthew G. Knepley 
725ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
7265cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
727ad0c61feSMatthew G. Knepley   PetscValidPointer(parent, 2);
728ad0c61feSMatthew G. Knepley   PetscValidPointer(name, 3);
729ad0c61feSMatthew G. Knepley   PetscValidPointer(value, 4);
730ad0c61feSMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
731ad0c61feSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
73260568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
73360568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
734f0b58479SMatthew G. Knepley   PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
735f0b58479SMatthew G. Knepley   if (datatype == PETSC_STRING) {
736f0b58479SMatthew G. Knepley     size_t len;
737f0b58479SMatthew G. Knepley 
738f0b58479SMatthew G. Knepley     PetscStackCallHDF5Return(len,H5Tget_size,(atype));
739f0b58479SMatthew G. Knepley     PetscStackCallHDF5(H5Tclose,(atype));
740f0b58479SMatthew G. Knepley     ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr);
741f0b58479SMatthew G. Knepley   }
74270efba33SBarry Smith   PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
743729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
74460568592SMatthew G. Knepley   PetscStackCallHDF5(H5Dclose,(obj));
745ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
746ad0c61feSMatthew G. Knepley }
747ad0c61feSMatthew G. Knepley 
7485cdeb108SMatthew G. Knepley static PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, const char name[], H5O_type_t otype, PetscBool *has)
7495cdeb108SMatthew G. Knepley {
7505cdeb108SMatthew G. Knepley   hid_t          h5;
7515cdeb108SMatthew G. Knepley   PetscErrorCode ierr;
7525cdeb108SMatthew G. Knepley 
7535cdeb108SMatthew G. Knepley   PetscFunctionBegin;
7545cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7555cdeb108SMatthew G. Knepley   PetscValidPointer(name, 2);
7565cdeb108SMatthew G. Knepley   PetscValidPointer(has, 3);
7575cdeb108SMatthew G. Knepley   *has = PETSC_FALSE;
758c0cd0301SJed Brown   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
7595cdeb108SMatthew G. Knepley   if (H5Lexists(h5, name, H5P_DEFAULT)) {
7605cdeb108SMatthew G. Knepley     H5O_info_t info;
7615cdeb108SMatthew G. Knepley     hid_t      obj;
7625cdeb108SMatthew G. Knepley 
763729ab6d9SBarry Smith     PetscStackCallHDF5Return(obj,H5Oopen,(h5, name, H5P_DEFAULT));
764729ab6d9SBarry Smith     PetscStackCallHDF5(H5Oget_info,(obj, &info));
7655cdeb108SMatthew G. Knepley     if (otype == info.type) *has = PETSC_TRUE;
766729ab6d9SBarry Smith     PetscStackCallHDF5(H5Oclose,(obj));
7675cdeb108SMatthew G. Knepley   }
7685cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
7695cdeb108SMatthew G. Knepley }
7705cdeb108SMatthew G. Knepley 
771df863907SAlex Fikl /*@C
772e7bdbf8eSMatthew G. Knepley  PetscViewerHDF5HasAttribute - Check whether a scalar attribute exists
773e7bdbf8eSMatthew G. Knepley 
774e7bdbf8eSMatthew G. Knepley   Input Parameters:
775e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer
776e7bdbf8eSMatthew G. Knepley . parent - The parent name
777e7bdbf8eSMatthew G. Knepley - name   - The attribute name
778e7bdbf8eSMatthew G. Knepley 
779e7bdbf8eSMatthew G. Knepley   Output Parameter:
780e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
781e7bdbf8eSMatthew G. Knepley 
782e7bdbf8eSMatthew G. Knepley   Level: advanced
783e7bdbf8eSMatthew G. Knepley 
784e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute()
785e7bdbf8eSMatthew G. Knepley @*/
786e7bdbf8eSMatthew G. Knepley PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
787e7bdbf8eSMatthew G. Knepley {
788729ab6d9SBarry Smith   hid_t          h5, dataset;
7895cdeb108SMatthew G. Knepley   htri_t         hhas;
7905cdeb108SMatthew G. Knepley   PetscBool      exists;
791e7bdbf8eSMatthew G. Knepley   PetscErrorCode ierr;
792e7bdbf8eSMatthew G. Knepley 
793e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
7945cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
795e7bdbf8eSMatthew G. Knepley   PetscValidPointer(parent, 2);
796e7bdbf8eSMatthew G. Knepley   PetscValidPointer(name, 3);
797e7bdbf8eSMatthew G. Knepley   PetscValidPointer(has, 4);
798e7bdbf8eSMatthew G. Knepley   *has = PETSC_FALSE;
799e7bdbf8eSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
8005cdeb108SMatthew G. Knepley   ierr = PetscViewerHDF5HasObject(viewer, parent, H5O_TYPE_DATASET, &exists);CHKERRQ(ierr);
8015cdeb108SMatthew G. Knepley   if (exists) {
802e7bdbf8eSMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
803729ab6d9SBarry Smith     PetscStackCall("H5Dopen2",dataset = H5Dopen2(h5, parent, H5P_DEFAULT));
804e7bdbf8eSMatthew G. Knepley #else
805729ab6d9SBarry Smith     PetscStackCall("H5Dopen",dataset = H5Dopen(h5, parent));
806e7bdbf8eSMatthew G. Knepley #endif
807729ab6d9SBarry Smith     if (dataset < 0) PetscFunctionReturn(0);
808729ab6d9SBarry Smith     PetscStackCall("H5Aexists",hhas = H5Aexists(dataset, name));
809729ab6d9SBarry Smith     if (hhas < 0) {
810729ab6d9SBarry Smith       PetscStackCallHDF5(H5Dclose,(dataset));
811729ab6d9SBarry Smith       PetscFunctionReturn(0);
812729ab6d9SBarry Smith     }
813729ab6d9SBarry Smith     PetscStackCallHDF5(H5Dclose,(dataset));
8145cdeb108SMatthew G. Knepley     *has = hhas ? PETSC_TRUE : PETSC_FALSE;
8155cdeb108SMatthew G. Knepley   }
816e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
817e7bdbf8eSMatthew G. Knepley }
818e7bdbf8eSMatthew G. Knepley 
81969a06e7bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadInitialize_Internal(PetscViewer viewer, const char name[], HDF5ReadCtx *ctx, PetscInt *timestep, PetscBool *complexVal)
820a5e1feadSVaclav Hapla {
82169a06e7bSVaclav Hapla   HDF5ReadCtx    h;
82269a06e7bSVaclav Hapla   PetscInt       timestep_;
82369a06e7bSVaclav Hapla   PetscBool      complexVal_ = PETSC_FALSE;
82469a06e7bSVaclav Hapla   const char    *groupname;
82569a06e7bSVaclav Hapla   char           vecgroup[PETSC_MAX_PATH_LEN];
826a5e1feadSVaclav Hapla   PetscErrorCode ierr;
827a5e1feadSVaclav Hapla 
828a5e1feadSVaclav Hapla   PetscFunctionBegin;
82969a06e7bSVaclav Hapla   ierr = PetscNew(&h);CHKERRQ(ierr);
83069a06e7bSVaclav Hapla   ierr = PetscViewerHDF5OpenGroup(viewer, &h->file, &h->group);CHKERRQ(ierr);
831a5e1feadSVaclav Hapla #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
83269a06e7bSVaclav Hapla   PetscStackCallHDF5Return(h->dataset,H5Dopen2,(h->group, name, H5P_DEFAULT));
833a5e1feadSVaclav Hapla #else
83469a06e7bSVaclav Hapla   PetscStackCallHDF5Return(h->dataset,H5Dopen,(h->group, name));
835a5e1feadSVaclav Hapla #endif
83669a06e7bSVaclav Hapla   PetscStackCallHDF5Return(h->dataspace,H5Dget_space,(h->dataset));
83769a06e7bSVaclav Hapla   ierr = PetscViewerHDF5GetTimestep(viewer, &timestep_);CHKERRQ(ierr);
83869a06e7bSVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer,&groupname);CHKERRQ(ierr);
83969a06e7bSVaclav Hapla   ierr = PetscSNPrintf(vecgroup,PETSC_MAX_PATH_LEN,"%s/%s",groupname,name);CHKERRQ(ierr);
84069a06e7bSVaclav Hapla   ierr = PetscViewerHDF5HasAttribute(viewer,vecgroup,"complex",&complexVal_);CHKERRQ(ierr);
8418af7cdc8SVaclav Hapla #if !defined(PETSC_USE_COMPLEX)
8428af7cdc8SVaclav Hapla   if (complexVal_) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"file contains complex numbers but PETSc not configured for them. Configure with --with-scalar-type=complex.");
8438af7cdc8SVaclav Hapla #endif
84469a06e7bSVaclav Hapla   *ctx = h;
84569a06e7bSVaclav Hapla   *timestep = timestep_;
84669a06e7bSVaclav Hapla   *complexVal = complexVal_;
84769a06e7bSVaclav Hapla   PetscFunctionReturn(0);
84869a06e7bSVaclav Hapla }
84969a06e7bSVaclav Hapla 
85069a06e7bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadFinalize_Internal(PetscViewer viewer, HDF5ReadCtx *ctx)
85169a06e7bSVaclav Hapla {
85269a06e7bSVaclav Hapla   HDF5ReadCtx    h;
85369a06e7bSVaclav Hapla   PetscErrorCode ierr;
85469a06e7bSVaclav Hapla 
85569a06e7bSVaclav Hapla   PetscFunctionBegin;
85669a06e7bSVaclav Hapla   h = *ctx;
85769a06e7bSVaclav Hapla   if (h->group != h->file) PetscStackCallHDF5(H5Gclose,(h->group));
85869a06e7bSVaclav Hapla   PetscStackCallHDF5(H5Sclose,(h->dataspace));
85969a06e7bSVaclav Hapla   PetscStackCallHDF5(H5Dclose,(h->dataset));
86069a06e7bSVaclav Hapla   ierr = PetscFree(*ctx);CHKERRQ(ierr);
86169a06e7bSVaclav Hapla   PetscFunctionReturn(0);
86269a06e7bSVaclav Hapla }
86369a06e7bSVaclav Hapla 
8648374c777SVaclav Hapla PetscErrorCode PetscViewerHDF5ReadSizes_Internal(PetscViewer viewer, HDF5ReadCtx ctx, PetscInt timestep, PetscBool complexVal, PetscLayout *map_)
86569a06e7bSVaclav Hapla {
86669a06e7bSVaclav Hapla   int            rdim, dim;
86769a06e7bSVaclav Hapla   hsize_t        dims[4];
868*45e21b7fSVaclav Hapla   PetscInt       bsInd, lenInd, bs, N;
869*45e21b7fSVaclav Hapla   PetscBool      dim2;
8708374c777SVaclav Hapla   PetscLayout    map;
8718374c777SVaclav Hapla   PetscErrorCode ierr;
87269a06e7bSVaclav Hapla 
87369a06e7bSVaclav Hapla   PetscFunctionBegin;
8748374c777SVaclav Hapla   if (!(*map_)) {
8758374c777SVaclav Hapla     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)viewer),map_);CHKERRQ(ierr);
8768374c777SVaclav Hapla   }
8778374c777SVaclav Hapla   map = *map_;
8789787f08bSVaclav Hapla   /* calculate expected number of dimensions */
879a5e1feadSVaclav Hapla   dim = 0;
880a5e1feadSVaclav Hapla   if (timestep >= 0) ++dim;
881a5e1feadSVaclav Hapla   ++dim; /* length in blocks */
882a5e1feadSVaclav Hapla   if (complexVal) ++dim;
8839787f08bSVaclav Hapla   /* get actual number of dimensions in dataset */
88469a06e7bSVaclav Hapla   PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(ctx->dataspace, dims, NULL));
8859787f08bSVaclav Hapla   /* calculate expected dimension indices */
8869787f08bSVaclav Hapla   lenInd = 0;
8879787f08bSVaclav Hapla   if (timestep >= 0) ++lenInd;
8889787f08bSVaclav Hapla   bsInd = lenInd + 1;
889*45e21b7fSVaclav Hapla   dim2 = PETSC_FALSE;
8909787f08bSVaclav Hapla   if (rdim == dim) {
891*45e21b7fSVaclav Hapla     bs = 1; /* support vectors stored as 1D array */
8929787f08bSVaclav Hapla   } else if (rdim == dim+1) {
893*45e21b7fSVaclav Hapla     bs = (PetscInt) dims[bsInd];
894*45e21b7fSVaclav Hapla     if (bs == 1) dim2 = PETSC_TRUE; /* vector with blocksize of 1, still stored as 2D array */
895a5e1feadSVaclav Hapla   } else {
8969787f08bSVaclav Hapla     SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Dimension of array in file %d not %d as expected", rdim, dim);
897a5e1feadSVaclav Hapla   }
898*45e21b7fSVaclav Hapla   N = (PetscInt) dims[lenInd]*bs;
8998374c777SVaclav Hapla 
9008374c777SVaclav Hapla   /* Set Vec sizes,blocksize,and type if not already set */
9018374c777SVaclav Hapla   if (map->bs < 0) {
902*45e21b7fSVaclav Hapla     ierr = PetscLayoutSetBlockSize(map, bs);CHKERRQ(ierr);
903*45e21b7fSVaclav 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);
9048374c777SVaclav Hapla   if (map->N < 0) {
905*45e21b7fSVaclav Hapla     ierr = PetscLayoutSetSize(map, N);CHKERRQ(ierr);
906*45e21b7fSVaclav 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);
90769a06e7bSVaclav Hapla   PetscFunctionReturn(0);
90869a06e7bSVaclav Hapla }
90969a06e7bSVaclav Hapla 
91069a06e7bSVaclav Hapla /* TODO DOC */
91169a06e7bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadSizes(PetscViewer viewer, const char name[], PetscInt *bs, PetscInt *N)
91269a06e7bSVaclav Hapla {
91369a06e7bSVaclav Hapla   HDF5ReadCtx    h;
91469a06e7bSVaclav Hapla   PetscInt       timestep;
91569a06e7bSVaclav Hapla   PetscBool      complexVal;
9168374c777SVaclav Hapla   PetscLayout    map=NULL;
91769a06e7bSVaclav Hapla   PetscErrorCode ierr;
91869a06e7bSVaclav Hapla 
91969a06e7bSVaclav Hapla   PetscFunctionBegin;
92069a06e7bSVaclav Hapla   ierr = PetscViewerHDF5ReadInitialize_Internal(viewer, name, &h, &timestep, &complexVal);CHKERRQ(ierr);
9218374c777SVaclav Hapla   ierr = PetscViewerHDF5ReadSizes_Internal(viewer, h, timestep, complexVal, &map);CHKERRQ(ierr);
92269a06e7bSVaclav Hapla   ierr = PetscViewerHDF5ReadFinalize_Internal(viewer, &h);CHKERRQ(ierr);
9238374c777SVaclav Hapla   if (bs) *bs = map->bs;
9248374c777SVaclav Hapla   if (N) *N = map->N;
9258374c777SVaclav Hapla   ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr);
926a5e1feadSVaclav Hapla   PetscFunctionReturn(0);
927a5e1feadSVaclav Hapla }
928a5e1feadSVaclav Hapla 
929a75e6a4aSMatthew G. Knepley /*
930a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
931a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
932a75e6a4aSMatthew G. Knepley */
933d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
934a75e6a4aSMatthew G. Knepley 
935a75e6a4aSMatthew G. Knepley /*@C
936a75e6a4aSMatthew G. Knepley   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
937a75e6a4aSMatthew G. Knepley 
938a75e6a4aSMatthew G. Knepley   Collective on MPI_Comm
939a75e6a4aSMatthew G. Knepley 
940a75e6a4aSMatthew G. Knepley   Input Parameter:
941a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer
942a75e6a4aSMatthew G. Knepley 
943a75e6a4aSMatthew G. Knepley   Level: intermediate
944a75e6a4aSMatthew G. Knepley 
945a75e6a4aSMatthew G. Knepley   Options Database Keys:
946a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name>
947a75e6a4aSMatthew G. Knepley 
948a75e6a4aSMatthew G. Knepley   Environmental variables:
949a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME
950a75e6a4aSMatthew G. Knepley 
951a75e6a4aSMatthew G. Knepley   Notes:
952a75e6a4aSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
953a75e6a4aSMatthew G. Knepley   an error code.  The HDF5 PetscViewer is usually used in the form
954a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
955a75e6a4aSMatthew G. Knepley 
956a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
957a75e6a4aSMatthew G. Knepley @*/
958a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
959a75e6a4aSMatthew G. Knepley {
960a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
961a75e6a4aSMatthew G. Knepley   PetscBool      flg;
962a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
963a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
964a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
965a75e6a4aSMatthew G. Knepley 
966a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
967a75e6a4aSMatthew 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);}
968a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
96912801b39SBarry Smith     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
970a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
971a75e6a4aSMatthew G. Knepley   }
97247435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
973a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
974a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
975a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
976a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
977a75e6a4aSMatthew G. Knepley     if (!flg) {
978a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname,"output.h5");
979a75e6a4aSMatthew G. Knepley       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
980a75e6a4aSMatthew G. Knepley     }
981a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
982a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
983a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
984a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
98547435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
986a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
987a75e6a4aSMatthew G. Knepley   }
988a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
989a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
990a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
991a75e6a4aSMatthew G. Knepley }
992