xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision b723ab350d0baf03346dec3580dd2ced69abe84a)
1af0996ceSBarry Smith #include <petsc/private/viewerimpl.h>    /*I   "petscsys.h"   I*/
2d70abbfaSBarry Smith #include <petscviewerhdf5.h>    /*I   "petscviewerhdf5.h"   I*/
3ded06c3fSVaclav Hapla #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE < 10800)
4ded06c3fSVaclav Hapla #error "PETSc needs HDF5 version >= 1.8.0"
5ded06c3fSVaclav Hapla #endif
65c6c1daeSBarry Smith 
7bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool*, H5O_type_t*);
806db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool*);
906db490cSVaclav Hapla 
105c6c1daeSBarry Smith typedef struct GroupList {
115c6c1daeSBarry Smith   const char       *name;
125c6c1daeSBarry Smith   struct GroupList *next;
135c6c1daeSBarry Smith } GroupList;
145c6c1daeSBarry Smith 
155c6c1daeSBarry Smith typedef struct {
165c6c1daeSBarry Smith   char          *filename;
175c6c1daeSBarry Smith   PetscFileMode btype;
185c6c1daeSBarry Smith   hid_t         file_id;
195c6c1daeSBarry Smith   PetscInt      timestep;
205c6c1daeSBarry Smith   GroupList     *groups;
2182ea9b62SBarry Smith   PetscBool     basedimension2;  /* save vectors and DMDA vectors with a dimension of at least 2 even if the bs/dof is 1 */
229a0502c6SHåkon Strandenes   PetscBool     spoutput;  /* write data in single precision even if PETSc is compiled with double precision PetscReal */
23058bd781SVaclav Hapla   char          *mataij_iname;
24058bd781SVaclav Hapla   char          *mataij_jname;
25058bd781SVaclav Hapla   char          *mataij_aname;
26058bd781SVaclav Hapla   char          *mataij_cname;
27*b723ab35SVaclav Hapla   PetscBool     mataij_names_set;
285c6c1daeSBarry Smith } PetscViewer_HDF5;
295c6c1daeSBarry Smith 
30eb5a92b4SVaclav Hapla struct _n_HDF5ReadCtx {
31eb5a92b4SVaclav Hapla   hid_t file, group, dataset, dataspace, plist;
32eb5a92b4SVaclav Hapla   PetscInt timestep;
3309dabeb0SVaclav Hapla   PetscBool complexVal, dim2, horizontal;
34eb5a92b4SVaclav Hapla };
35eb5a92b4SVaclav Hapla typedef struct _n_HDF5ReadCtx* HDF5ReadCtx;
36eb5a92b4SVaclav Hapla 
376c132bc1SVaclav Hapla static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char objname[], char **fullpath)
386c132bc1SVaclav Hapla {
396c132bc1SVaclav Hapla   const char *group;
406c132bc1SVaclav Hapla   char buf[PETSC_MAX_PATH_LEN]="";
416c132bc1SVaclav Hapla   PetscErrorCode ierr;
426c132bc1SVaclav Hapla 
436c132bc1SVaclav Hapla   PetscFunctionBegin;
446c132bc1SVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr);
456c132bc1SVaclav Hapla   ierr = PetscStrcat(buf, group);CHKERRQ(ierr);
466c132bc1SVaclav Hapla   ierr = PetscStrcat(buf, "/");CHKERRQ(ierr);
476c132bc1SVaclav Hapla   ierr = PetscStrcat(buf, objname);CHKERRQ(ierr);
486c132bc1SVaclav Hapla   ierr = PetscStrallocpy(buf, fullpath);CHKERRQ(ierr);
496c132bc1SVaclav Hapla   PetscFunctionReturn(0);
506c132bc1SVaclav Hapla }
516c132bc1SVaclav Hapla 
52577e0e04SVaclav Hapla static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj)
53577e0e04SVaclav Hapla {
54577e0e04SVaclav Hapla   PetscBool has;
55577e0e04SVaclav Hapla   const char *group;
56577e0e04SVaclav Hapla   PetscErrorCode ierr;
57577e0e04SVaclav Hapla 
58577e0e04SVaclav Hapla   PetscFunctionBegin;
59577e0e04SVaclav Hapla   if (!obj->name) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
60577e0e04SVaclav Hapla   ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr);
61577e0e04SVaclav Hapla   if (!has) {
62577e0e04SVaclav Hapla     ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr);
63577e0e04SVaclav Hapla     SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) %s not stored in group %s", obj->name, group);
64577e0e04SVaclav Hapla   }
65577e0e04SVaclav Hapla   PetscFunctionReturn(0);
66577e0e04SVaclav Hapla }
67577e0e04SVaclav Hapla 
684416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v)
6982ea9b62SBarry Smith {
7082ea9b62SBarry Smith   PetscErrorCode   ierr;
7182ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;
7282ea9b62SBarry Smith 
7382ea9b62SBarry Smith   PetscFunctionBegin;
7482ea9b62SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr);
7582ea9b62SBarry Smith   ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr);
769a0502c6SHåkon Strandenes   ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr);
7782ea9b62SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
7882ea9b62SBarry Smith   PetscFunctionReturn(0);
7982ea9b62SBarry Smith }
8082ea9b62SBarry Smith 
815c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
825c6c1daeSBarry Smith {
835c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
845c6c1daeSBarry Smith   PetscErrorCode   ierr;
855c6c1daeSBarry Smith 
865c6c1daeSBarry Smith   PetscFunctionBegin;
875c6c1daeSBarry Smith   ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);
88729ab6d9SBarry Smith   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
895c6c1daeSBarry Smith   PetscFunctionReturn(0);
905c6c1daeSBarry Smith }
915c6c1daeSBarry Smith 
925c6c1daeSBarry Smith PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
935c6c1daeSBarry Smith {
945c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
955c6c1daeSBarry Smith   PetscErrorCode   ierr;
965c6c1daeSBarry Smith 
975c6c1daeSBarry Smith   PetscFunctionBegin;
985c6c1daeSBarry Smith   ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr);
995c6c1daeSBarry Smith   while (hdf5->groups) {
1005c6c1daeSBarry Smith     GroupList *tmp = hdf5->groups->next;
1015c6c1daeSBarry Smith 
1025c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups->name);CHKERRQ(ierr);
1035c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups);CHKERRQ(ierr);
1045c6c1daeSBarry Smith     hdf5->groups = tmp;
1055c6c1daeSBarry Smith   }
10661912a0cSVaclav Hapla   ierr = PetscFree(hdf5->mataij_iname);CHKERRQ(ierr);
10761912a0cSVaclav Hapla   ierr = PetscFree(hdf5->mataij_jname);CHKERRQ(ierr);
10861912a0cSVaclav Hapla   ierr = PetscFree(hdf5->mataij_aname);CHKERRQ(ierr);
10961912a0cSVaclav Hapla   ierr = PetscFree(hdf5->mataij_cname);CHKERRQ(ierr);
1105c6c1daeSBarry Smith   ierr = PetscFree(hdf5);CHKERRQ(ierr);
1110b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
112d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr);
1130b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
114058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr);
115058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr);
116058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetAIJNames_C",NULL);CHKERRQ(ierr);
117058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5GetAIJNames_C",NULL);CHKERRQ(ierr);
1185c6c1daeSBarry Smith   PetscFunctionReturn(0);
1195c6c1daeSBarry Smith }
1205c6c1daeSBarry Smith 
1215c6c1daeSBarry Smith PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
1225c6c1daeSBarry Smith {
1235c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1245c6c1daeSBarry Smith 
1255c6c1daeSBarry Smith   PetscFunctionBegin;
1265c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1275c6c1daeSBarry Smith   hdf5->btype = type;
1285c6c1daeSBarry Smith   PetscFunctionReturn(0);
1295c6c1daeSBarry Smith }
1305c6c1daeSBarry Smith 
13182ea9b62SBarry Smith PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
13282ea9b62SBarry Smith {
13382ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
13482ea9b62SBarry Smith 
13582ea9b62SBarry Smith   PetscFunctionBegin;
13682ea9b62SBarry Smith   hdf5->basedimension2 = flg;
13782ea9b62SBarry Smith   PetscFunctionReturn(0);
13882ea9b62SBarry Smith }
13982ea9b62SBarry Smith 
140df863907SAlex Fikl /*@
14182ea9b62SBarry Smith      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
14282ea9b62SBarry Smith        dimension of 2.
14382ea9b62SBarry Smith 
14482ea9b62SBarry Smith     Logically Collective on PetscViewer
14582ea9b62SBarry Smith 
14682ea9b62SBarry Smith   Input Parameters:
14782ea9b62SBarry Smith +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
14882ea9b62SBarry 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
14982ea9b62SBarry Smith 
15082ea9b62SBarry Smith   Options Database:
15182ea9b62SBarry 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
15282ea9b62SBarry Smith 
15382ea9b62SBarry Smith 
15495452b02SPatrick Sanan   Notes:
15595452b02SPatrick 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
15682ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
15782ea9b62SBarry Smith 
15882ea9b62SBarry Smith   Level: intermediate
15982ea9b62SBarry Smith 
16082ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
16182ea9b62SBarry Smith 
16282ea9b62SBarry Smith @*/
16382ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
16482ea9b62SBarry Smith {
16582ea9b62SBarry Smith   PetscErrorCode ierr;
16682ea9b62SBarry Smith 
16782ea9b62SBarry Smith   PetscFunctionBegin;
16882ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
16982ea9b62SBarry Smith   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
17082ea9b62SBarry Smith   PetscFunctionReturn(0);
17182ea9b62SBarry Smith }
17282ea9b62SBarry Smith 
173df863907SAlex Fikl /*@
17482ea9b62SBarry Smith      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
17582ea9b62SBarry Smith        dimension of 2.
17682ea9b62SBarry Smith 
17782ea9b62SBarry Smith     Logically Collective on PetscViewer
17882ea9b62SBarry Smith 
17982ea9b62SBarry Smith   Input Parameter:
18082ea9b62SBarry Smith .  viewer - the PetscViewer, must be of type HDF5
18182ea9b62SBarry Smith 
18282ea9b62SBarry Smith   Output Parameter:
18382ea9b62SBarry 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
18482ea9b62SBarry Smith 
18595452b02SPatrick Sanan   Notes:
18695452b02SPatrick 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
18782ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
18882ea9b62SBarry Smith 
18982ea9b62SBarry Smith   Level: intermediate
19082ea9b62SBarry Smith 
19182ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
19282ea9b62SBarry Smith 
19382ea9b62SBarry Smith @*/
19482ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
19582ea9b62SBarry Smith {
19682ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
19782ea9b62SBarry Smith 
19882ea9b62SBarry Smith   PetscFunctionBegin;
19982ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
20082ea9b62SBarry Smith   *flg = hdf5->basedimension2;
20182ea9b62SBarry Smith   PetscFunctionReturn(0);
20282ea9b62SBarry Smith }
20382ea9b62SBarry Smith 
2049a0502c6SHåkon Strandenes PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
2059a0502c6SHåkon Strandenes {
2069a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2079a0502c6SHåkon Strandenes 
2089a0502c6SHåkon Strandenes   PetscFunctionBegin;
2099a0502c6SHåkon Strandenes   hdf5->spoutput = flg;
2109a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2119a0502c6SHåkon Strandenes }
2129a0502c6SHåkon Strandenes 
213df863907SAlex Fikl /*@
2149a0502c6SHåkon Strandenes      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
2159a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2169a0502c6SHåkon Strandenes 
2179a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2189a0502c6SHåkon Strandenes 
2199a0502c6SHåkon Strandenes   Input Parameters:
2209a0502c6SHåkon Strandenes +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
2219a0502c6SHåkon Strandenes -  flg - if PETSC_TRUE the data will be written to disk with single precision
2229a0502c6SHåkon Strandenes 
2239a0502c6SHåkon Strandenes   Options Database:
2249a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
2259a0502c6SHåkon Strandenes 
2269a0502c6SHåkon Strandenes 
22795452b02SPatrick Sanan   Notes:
22895452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2299a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2309a0502c6SHåkon Strandenes 
2319a0502c6SHåkon Strandenes   Level: intermediate
2329a0502c6SHåkon Strandenes 
2339a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2349a0502c6SHåkon Strandenes           PetscReal
2359a0502c6SHåkon Strandenes 
2369a0502c6SHåkon Strandenes @*/
2379a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
2389a0502c6SHåkon Strandenes {
2399a0502c6SHåkon Strandenes   PetscErrorCode ierr;
2409a0502c6SHåkon Strandenes 
2419a0502c6SHåkon Strandenes   PetscFunctionBegin;
2429a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2439a0502c6SHåkon Strandenes   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
2449a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2459a0502c6SHåkon Strandenes }
2469a0502c6SHåkon Strandenes 
247df863907SAlex Fikl /*@
2489a0502c6SHåkon Strandenes      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
2499a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2509a0502c6SHåkon Strandenes 
2519a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2529a0502c6SHåkon Strandenes 
2539a0502c6SHåkon Strandenes   Input Parameter:
2549a0502c6SHåkon Strandenes .  viewer - the PetscViewer, must be of type HDF5
2559a0502c6SHåkon Strandenes 
2569a0502c6SHåkon Strandenes   Output Parameter:
2579a0502c6SHåkon Strandenes .  flg - if PETSC_TRUE the data will be written to disk with single precision
2589a0502c6SHåkon Strandenes 
25995452b02SPatrick Sanan   Notes:
26095452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2619a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2629a0502c6SHåkon Strandenes 
2639a0502c6SHåkon Strandenes   Level: intermediate
2649a0502c6SHåkon Strandenes 
2659a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2669a0502c6SHåkon Strandenes           PetscReal
2679a0502c6SHåkon Strandenes 
2689a0502c6SHåkon Strandenes @*/
2699a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
2709a0502c6SHåkon Strandenes {
2719a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2729a0502c6SHåkon Strandenes 
2739a0502c6SHåkon Strandenes   PetscFunctionBegin;
2749a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2759a0502c6SHåkon Strandenes   *flg = hdf5->spoutput;
2769a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2779a0502c6SHåkon Strandenes }
2789a0502c6SHåkon Strandenes 
2795c6c1daeSBarry Smith PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
2805c6c1daeSBarry Smith {
2815c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2825c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
2835c6c1daeSBarry Smith   MPI_Info          info = MPI_INFO_NULL;
2845c6c1daeSBarry Smith #endif
2855c6c1daeSBarry Smith   hid_t             plist_id;
2865c6c1daeSBarry Smith   PetscErrorCode    ierr;
2875c6c1daeSBarry Smith 
2885c6c1daeSBarry Smith   PetscFunctionBegin;
289f683cc58SToby Isaac   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
290f683cc58SToby Isaac   if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);}
2915c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
2925c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
293729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
2945c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
295729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info));
2965c6c1daeSBarry Smith #endif
2975c6c1daeSBarry Smith   /* Create or open the file collectively */
2985c6c1daeSBarry Smith   switch (hdf5->btype) {
2995c6c1daeSBarry Smith   case FILE_MODE_READ:
300729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
3015c6c1daeSBarry Smith     break;
3025c6c1daeSBarry Smith   case FILE_MODE_APPEND:
303729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
3045c6c1daeSBarry Smith     break;
3055c6c1daeSBarry Smith   case FILE_MODE_WRITE:
306729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
3075c6c1daeSBarry Smith     break;
3085c6c1daeSBarry Smith   default:
3095c6c1daeSBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
3105c6c1daeSBarry Smith   }
3115c6c1daeSBarry Smith   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
312729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
3135c6c1daeSBarry Smith   PetscFunctionReturn(0);
3145c6c1daeSBarry Smith }
3155c6c1daeSBarry Smith 
316d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
317d1232d7fSToby Isaac {
318d1232d7fSToby Isaac   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
319d1232d7fSToby Isaac 
320d1232d7fSToby Isaac   PetscFunctionBegin;
321d1232d7fSToby Isaac   *name = vhdf5->filename;
322d1232d7fSToby Isaac   PetscFunctionReturn(0);
323d1232d7fSToby Isaac }
324d1232d7fSToby Isaac 
325058bd781SVaclav Hapla PetscErrorCode  PetscViewerHDF5SetAIJNames_HDF5(PetscViewer viewer, const char iname[], const char jname[], const char aname[], const char cname[])
326058bd781SVaclav Hapla {
327058bd781SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
328058bd781SVaclav Hapla   PetscErrorCode ierr;
329058bd781SVaclav Hapla 
330058bd781SVaclav Hapla   PetscFunctionBegin;
331058bd781SVaclav Hapla   ierr = PetscFree(hdf5->mataij_iname);CHKERRQ(ierr);
332058bd781SVaclav Hapla   ierr = PetscFree(hdf5->mataij_jname);CHKERRQ(ierr);
333058bd781SVaclav Hapla   ierr = PetscFree(hdf5->mataij_aname);CHKERRQ(ierr);
334058bd781SVaclav Hapla   ierr = PetscFree(hdf5->mataij_cname);CHKERRQ(ierr);
335058bd781SVaclav Hapla   ierr = PetscStrallocpy(iname,&hdf5->mataij_iname);CHKERRQ(ierr);
336058bd781SVaclav Hapla   ierr = PetscStrallocpy(jname,&hdf5->mataij_jname);CHKERRQ(ierr);
337058bd781SVaclav Hapla   ierr = PetscStrallocpy(aname,&hdf5->mataij_aname);CHKERRQ(ierr);
338058bd781SVaclav Hapla   ierr = PetscStrallocpy(cname,&hdf5->mataij_cname);CHKERRQ(ierr);
339*b723ab35SVaclav Hapla   hdf5->mataij_names_set = PETSC_TRUE;
340058bd781SVaclav Hapla   PetscFunctionReturn(0);
341058bd781SVaclav Hapla }
342058bd781SVaclav Hapla 
343058bd781SVaclav Hapla /*@C
344058bd781SVaclav Hapla   PetscViewerHDF5SetAIJNames - Set the names of the datasets representing the three AIJ (CRS) arrays and the name of the attribute storing the number of columns within the HDF5 file.
345058bd781SVaclav Hapla 
346058bd781SVaclav Hapla   Collective on PetscViewer
347058bd781SVaclav Hapla 
348058bd781SVaclav Hapla   Input Parameters:
349058bd781SVaclav Hapla +  viewer - the PetscViewer; either ASCII or binary
350058bd781SVaclav Hapla .  iname - name of dataset i representing row pointers; that is i[0] = 0, i[row] = i[row-1] + number of elements in that row of the matrix
351058bd781SVaclav Hapla .  jname - name of dataset j representing column indices
352058bd781SVaclav Hapla .  aname - name of dataset a representing matrix values
353058bd781SVaclav Hapla -  cname - name of attribute stoting column count
354058bd781SVaclav Hapla 
355058bd781SVaclav Hapla   Level: advanced
356058bd781SVaclav Hapla 
357058bd781SVaclav Hapla   Notes:
358058bd781SVaclav Hapla   Current defaults are (iname, jname, aname, cname) = ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be readily loaded.
359058bd781SVaclav Hapla 
360058bd781SVaclav Hapla .seealso: MatLoad(), PetscViewerCreate(), PetscViewerSetType(), PETSCVIEWERHDF5, PetscViewerHDF5GetAIJNames()
361058bd781SVaclav Hapla @*/
362058bd781SVaclav Hapla PetscErrorCode  PetscViewerHDF5SetAIJNames(PetscViewer viewer, const char iname[], const char jname[], const char aname[], const char cname[])
363058bd781SVaclav Hapla {
364058bd781SVaclav Hapla   PetscErrorCode ierr;
365058bd781SVaclav Hapla 
366058bd781SVaclav Hapla   PetscFunctionBegin;
367058bd781SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
368058bd781SVaclav Hapla   PetscValidCharPointer(iname,2);
369058bd781SVaclav Hapla   PetscValidCharPointer(jname,3);
370058bd781SVaclav Hapla   PetscValidCharPointer(aname,4);
371058bd781SVaclav Hapla   PetscValidCharPointer(cname,5);
372058bd781SVaclav Hapla   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetAIJNames_C",(PetscViewer,const char[],const char[],const char[],const char[]),(viewer,iname,jname,aname,cname));CHKERRQ(ierr);
373058bd781SVaclav Hapla   PetscFunctionReturn(0);
374058bd781SVaclav Hapla }
375058bd781SVaclav Hapla 
376058bd781SVaclav Hapla PetscErrorCode  PetscViewerHDF5GetAIJNames_HDF5(PetscViewer viewer, const char *iname[], const char *jname[], const char *aname[], const char *cname[])
377058bd781SVaclav Hapla {
378058bd781SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
379058bd781SVaclav Hapla 
380058bd781SVaclav Hapla   PetscFunctionBegin;
381058bd781SVaclav Hapla   *iname = hdf5->mataij_iname;
382058bd781SVaclav Hapla   *jname = hdf5->mataij_jname;
383058bd781SVaclav Hapla   *aname = hdf5->mataij_aname;
384058bd781SVaclav Hapla   *cname = hdf5->mataij_cname;
385058bd781SVaclav Hapla   PetscFunctionReturn(0);
386058bd781SVaclav Hapla }
387058bd781SVaclav Hapla 
388058bd781SVaclav Hapla /*@C
389058bd781SVaclav Hapla   PetscViewerHDF5GetAIJNames - Get the names of the datasets representing the three AIJ (CRS) arrays and the name of the attribute storing the number of columns within the HDF5 file.
390058bd781SVaclav Hapla 
391058bd781SVaclav Hapla   Collective on PetscViewer
392058bd781SVaclav Hapla 
393058bd781SVaclav Hapla   Input Parameters:
394058bd781SVaclav Hapla .  viewer - the PetscViewer; either ASCII or binary
395058bd781SVaclav Hapla 
396058bd781SVaclav Hapla   Output Parameters:
397058bd781SVaclav Hapla +  iname - name of dataset i representing row pointers; that is i[0] = 0, i[row] = i[row-1] + number of elements in that row of the matrix
398058bd781SVaclav Hapla .  jname - name of dataset j representing column indices
399058bd781SVaclav Hapla .  aname - name of dataset a representing matrix values
400058bd781SVaclav Hapla -  cname - name of attribute stoting column count
401058bd781SVaclav Hapla 
402058bd781SVaclav Hapla   Level: advanced
403058bd781SVaclav Hapla 
404058bd781SVaclav Hapla   Notes:
405058bd781SVaclav Hapla   Current defaults are (iname, jname, aname, cname) = ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be readily loaded.
406058bd781SVaclav Hapla 
407058bd781SVaclav Hapla .seealso: MatLoad(), PetscViewerCreate(), PetscViewerSetType(), PETSCVIEWERHDF5, PetscViewerHDF5SetAIJNames()
408058bd781SVaclav Hapla @*/
409058bd781SVaclav Hapla PetscErrorCode  PetscViewerHDF5GetAIJNames(PetscViewer viewer, const char *iname[], const char *jname[], const char *aname[], const char *cname[])
410058bd781SVaclav Hapla {
411058bd781SVaclav Hapla   PetscErrorCode ierr;
412058bd781SVaclav Hapla 
413058bd781SVaclav Hapla   PetscFunctionBegin;
414058bd781SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
415058bd781SVaclav Hapla   PetscValidPointer(iname,2);
416058bd781SVaclav Hapla   PetscValidPointer(jname,3);
417058bd781SVaclav Hapla   PetscValidPointer(aname,4);
418058bd781SVaclav Hapla   PetscValidPointer(cname,5);
419058bd781SVaclav Hapla   ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetAIJNames_C",(PetscViewer,const char*[],const char*[],const char*[],const char*[]),(viewer,iname,jname,aname,cname));CHKERRQ(ierr);
420058bd781SVaclav Hapla   PetscFunctionReturn(0);
421058bd781SVaclav Hapla }
422058bd781SVaclav Hapla 
423*b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
424*b723ab35SVaclav Hapla {
425*b723ab35SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
426*b723ab35SVaclav Hapla   PetscErrorCode   ierr;
427*b723ab35SVaclav Hapla 
428*b723ab35SVaclav Hapla   PetscFunctionBegin;
429*b723ab35SVaclav Hapla   if (!hdf5->mataij_names_set) {
430*b723ab35SVaclav Hapla     ierr = PetscViewerHDF5SetAIJNames_HDF5(viewer,"jc","ir","data","MATLAB_sparse");CHKERRQ(ierr);
431*b723ab35SVaclav Hapla   }
432*b723ab35SVaclav Hapla   PetscFunctionReturn(0);
433*b723ab35SVaclav Hapla }
434*b723ab35SVaclav Hapla 
4358556b5ebSBarry Smith /*MC
4368556b5ebSBarry Smith    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
4378556b5ebSBarry Smith 
4388556b5ebSBarry Smith 
4398556b5ebSBarry Smith .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
4408556b5ebSBarry Smith            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
4418556b5ebSBarry Smith            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
4428556b5ebSBarry Smith            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
4438556b5ebSBarry Smith 
4441b266c99SBarry Smith   Level: beginner
4458556b5ebSBarry Smith M*/
446d1232d7fSToby Isaac 
4478cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
4485c6c1daeSBarry Smith {
4495c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
4505c6c1daeSBarry Smith   PetscErrorCode   ierr;
4515c6c1daeSBarry Smith 
4525c6c1daeSBarry Smith   PetscFunctionBegin;
453b00a9115SJed Brown   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
4545c6c1daeSBarry Smith 
4555c6c1daeSBarry Smith   v->data                = (void*) hdf5;
4565c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
45782ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
458*b723ab35SVaclav Hapla   v->ops->setup          = PetscViewerSetUp_HDF5;
4595c6c1daeSBarry Smith   v->ops->flush          = 0;
4605c6c1daeSBarry Smith   hdf5->btype            = (PetscFileMode) -1;
4615c6c1daeSBarry Smith   hdf5->filename         = 0;
4625c6c1daeSBarry Smith   hdf5->timestep         = -1;
4630298fd71SBarry Smith   hdf5->groups           = NULL;
4645c6c1daeSBarry Smith 
465*b723ab35SVaclav Hapla   hdf5->mataij_iname     = NULL;
466*b723ab35SVaclav Hapla   hdf5->mataij_jname     = NULL;
467*b723ab35SVaclav Hapla   hdf5->mataij_aname     = NULL;
468*b723ab35SVaclav Hapla   hdf5->mataij_cname     = NULL;
469*b723ab35SVaclav Hapla   hdf5->mataij_names_set = PETSC_FALSE;
470058bd781SVaclav Hapla 
4710b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
472d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr);
4730b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
47482ea9b62SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
4759a0502c6SHåkon Strandenes   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr);
476058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetAIJNames_C",PetscViewerHDF5SetAIJNames_HDF5);CHKERRQ(ierr);
477058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetAIJNames_C",PetscViewerHDF5GetAIJNames_HDF5);CHKERRQ(ierr);
4785c6c1daeSBarry Smith   PetscFunctionReturn(0);
4795c6c1daeSBarry Smith }
4805c6c1daeSBarry Smith 
4815c6c1daeSBarry Smith /*@C
4825c6c1daeSBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
4835c6c1daeSBarry Smith 
4845c6c1daeSBarry Smith    Collective on MPI_Comm
4855c6c1daeSBarry Smith 
4865c6c1daeSBarry Smith    Input Parameters:
4875c6c1daeSBarry Smith +  comm - MPI communicator
4885c6c1daeSBarry Smith .  name - name of file
4895c6c1daeSBarry Smith -  type - type of file
4905c6c1daeSBarry Smith $    FILE_MODE_WRITE - create new file for binary output
4915c6c1daeSBarry Smith $    FILE_MODE_READ - open existing file for binary input
4925c6c1daeSBarry Smith $    FILE_MODE_APPEND - open existing file for binary output
4935c6c1daeSBarry Smith 
4945c6c1daeSBarry Smith    Output Parameter:
4955c6c1daeSBarry Smith .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
4965c6c1daeSBarry Smith 
49782ea9b62SBarry Smith   Options Database:
49882ea9b62SBarry 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
4999a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
50082ea9b62SBarry Smith 
5015c6c1daeSBarry Smith    Level: beginner
5025c6c1daeSBarry Smith 
5035c6c1daeSBarry Smith    Note:
5045c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
5055c6c1daeSBarry Smith 
5065c6c1daeSBarry Smith    Concepts: HDF5 files
5075c6c1daeSBarry Smith    Concepts: PetscViewerHDF5^creating
5085c6c1daeSBarry Smith 
5096a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
5109a0502c6SHåkon Strandenes           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
511a56f64adSBarry Smith           MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
5125c6c1daeSBarry Smith @*/
5135c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
5145c6c1daeSBarry Smith {
5155c6c1daeSBarry Smith   PetscErrorCode ierr;
5165c6c1daeSBarry Smith 
5175c6c1daeSBarry Smith   PetscFunctionBegin;
5185c6c1daeSBarry Smith   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
5195c6c1daeSBarry Smith   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
5205c6c1daeSBarry Smith   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
5215c6c1daeSBarry Smith   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
5225c6c1daeSBarry Smith   PetscFunctionReturn(0);
5235c6c1daeSBarry Smith }
5245c6c1daeSBarry Smith 
5255c6c1daeSBarry Smith /*@C
5265c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
5275c6c1daeSBarry Smith 
5285c6c1daeSBarry Smith   Not collective
5295c6c1daeSBarry Smith 
5305c6c1daeSBarry Smith   Input Parameter:
5315c6c1daeSBarry Smith . viewer - the PetscViewer
5325c6c1daeSBarry Smith 
5335c6c1daeSBarry Smith   Output Parameter:
5345c6c1daeSBarry Smith . file_id - The file id
5355c6c1daeSBarry Smith 
5365c6c1daeSBarry Smith   Level: intermediate
5375c6c1daeSBarry Smith 
5385c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open()
5395c6c1daeSBarry Smith @*/
5405c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
5415c6c1daeSBarry Smith {
5425c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5435c6c1daeSBarry Smith 
5445c6c1daeSBarry Smith   PetscFunctionBegin;
5455c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5465c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
5475c6c1daeSBarry Smith   PetscFunctionReturn(0);
5485c6c1daeSBarry Smith }
5495c6c1daeSBarry Smith 
5505c6c1daeSBarry Smith /*@C
5515c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
5525c6c1daeSBarry Smith 
5535c6c1daeSBarry Smith   Not collective
5545c6c1daeSBarry Smith 
5555c6c1daeSBarry Smith   Input Parameters:
5565c6c1daeSBarry Smith + viewer - the PetscViewer
5575c6c1daeSBarry Smith - name - The group name
5585c6c1daeSBarry Smith 
5595c6c1daeSBarry Smith   Level: intermediate
5605c6c1daeSBarry Smith 
561874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
5625c6c1daeSBarry Smith @*/
5635c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
5645c6c1daeSBarry Smith {
5655c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5665c6c1daeSBarry Smith   GroupList        *groupNode;
5675c6c1daeSBarry Smith   PetscErrorCode   ierr;
5685c6c1daeSBarry Smith 
5695c6c1daeSBarry Smith   PetscFunctionBegin;
5705c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5715c6c1daeSBarry Smith   PetscValidCharPointer(name,2);
57295dccacaSBarry Smith   ierr = PetscNew(&groupNode);CHKERRQ(ierr);
5735c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
574a297a907SKarl Rupp 
5755c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
5765c6c1daeSBarry Smith   hdf5->groups    = groupNode;
5775c6c1daeSBarry Smith   PetscFunctionReturn(0);
5785c6c1daeSBarry Smith }
5795c6c1daeSBarry Smith 
5803ef9c667SSatish Balay /*@
5815c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
5825c6c1daeSBarry Smith 
5835c6c1daeSBarry Smith   Not collective
5845c6c1daeSBarry Smith 
5855c6c1daeSBarry Smith   Input Parameter:
5865c6c1daeSBarry Smith . viewer - the PetscViewer
5875c6c1daeSBarry Smith 
5885c6c1daeSBarry Smith   Level: intermediate
5895c6c1daeSBarry Smith 
590874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
5915c6c1daeSBarry Smith @*/
5925c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
5935c6c1daeSBarry Smith {
5945c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5955c6c1daeSBarry Smith   GroupList        *groupNode;
5965c6c1daeSBarry Smith   PetscErrorCode   ierr;
5975c6c1daeSBarry Smith 
5985c6c1daeSBarry Smith   PetscFunctionBegin;
5995c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
60082f516ccSBarry Smith   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
6015c6c1daeSBarry Smith   groupNode    = hdf5->groups;
6025c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
6035c6c1daeSBarry Smith   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
6045c6c1daeSBarry Smith   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
6055c6c1daeSBarry Smith   PetscFunctionReturn(0);
6065c6c1daeSBarry Smith }
6075c6c1daeSBarry Smith 
6085c6c1daeSBarry Smith /*@C
609874270d9SVaclav Hapla   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup().
610874270d9SVaclav Hapla   If none has been assigned, returns NULL.
6115c6c1daeSBarry Smith 
6125c6c1daeSBarry Smith   Not collective
6135c6c1daeSBarry Smith 
6145c6c1daeSBarry Smith   Input Parameter:
6155c6c1daeSBarry Smith . viewer - the PetscViewer
6165c6c1daeSBarry Smith 
6175c6c1daeSBarry Smith   Output Parameter:
6185c6c1daeSBarry Smith . name - The group name
6195c6c1daeSBarry Smith 
6205c6c1daeSBarry Smith   Level: intermediate
6215c6c1daeSBarry Smith 
622874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
6235c6c1daeSBarry Smith @*/
6245c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
6255c6c1daeSBarry Smith {
6265c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
6275c6c1daeSBarry Smith 
6285c6c1daeSBarry Smith   PetscFunctionBegin;
6295c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
630c959eef4SJed Brown   PetscValidPointer(name,2);
631a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
6320298fd71SBarry Smith   else *name = NULL;
6335c6c1daeSBarry Smith   PetscFunctionReturn(0);
6345c6c1daeSBarry Smith }
6355c6c1daeSBarry Smith 
6368c557b5aSVaclav Hapla /*@
637874270d9SVaclav Hapla   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(),
638874270d9SVaclav Hapla   and return this group's ID and file ID.
639874270d9SVaclav Hapla   If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID.
640874270d9SVaclav Hapla 
641874270d9SVaclav Hapla   Not collective
642874270d9SVaclav Hapla 
643874270d9SVaclav Hapla   Input Parameter:
644874270d9SVaclav Hapla . viewer - the PetscViewer
645874270d9SVaclav Hapla 
646874270d9SVaclav Hapla   Output Parameter:
647874270d9SVaclav Hapla + fileId - The HDF5 file ID
648874270d9SVaclav Hapla - groupId - The HDF5 group ID
649874270d9SVaclav Hapla 
650874270d9SVaclav Hapla   Level: intermediate
651874270d9SVaclav Hapla 
652874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
653874270d9SVaclav Hapla @*/
65454dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
65554dbf706SVaclav Hapla {
65690e03537SVaclav Hapla   hid_t          file_id;
65790e03537SVaclav Hapla   H5O_type_t     type;
65854dbf706SVaclav Hapla   const char     *groupName = NULL;
65954dbf706SVaclav Hapla   PetscErrorCode ierr;
66054dbf706SVaclav Hapla 
66154dbf706SVaclav Hapla   PetscFunctionBegin;
66254dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr);
66354dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr);
66490e03537SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, PETSC_TRUE, NULL, &type);CHKERRQ(ierr);
66590e03537SVaclav Hapla   if (type != H5O_TYPE_GROUP) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s resolves to something which is not a group", groupName);
66690e03537SVaclav Hapla   PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT));
66754dbf706SVaclav Hapla   *fileId  = file_id;
66854dbf706SVaclav Hapla   PetscFunctionReturn(0);
66954dbf706SVaclav Hapla }
67054dbf706SVaclav Hapla 
6713ef9c667SSatish Balay /*@
6725c6c1daeSBarry Smith   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
6735c6c1daeSBarry Smith 
6745c6c1daeSBarry Smith   Not collective
6755c6c1daeSBarry Smith 
6765c6c1daeSBarry Smith   Input Parameter:
6775c6c1daeSBarry Smith . viewer - the PetscViewer
6785c6c1daeSBarry Smith 
6795c6c1daeSBarry Smith   Level: intermediate
6805c6c1daeSBarry Smith 
6815c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
6825c6c1daeSBarry Smith @*/
6835c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
6845c6c1daeSBarry Smith {
6855c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
6865c6c1daeSBarry Smith 
6875c6c1daeSBarry Smith   PetscFunctionBegin;
6885c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
6895c6c1daeSBarry Smith   ++hdf5->timestep;
6905c6c1daeSBarry Smith   PetscFunctionReturn(0);
6915c6c1daeSBarry Smith }
6925c6c1daeSBarry Smith 
6933ef9c667SSatish Balay /*@
6945c6c1daeSBarry Smith   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
6955c6c1daeSBarry Smith   of -1 disables blocking with timesteps.
6965c6c1daeSBarry Smith 
6975c6c1daeSBarry Smith   Not collective
6985c6c1daeSBarry Smith 
6995c6c1daeSBarry Smith   Input Parameters:
7005c6c1daeSBarry Smith + viewer - the PetscViewer
7015c6c1daeSBarry Smith - timestep - The timestep number
7025c6c1daeSBarry Smith 
7035c6c1daeSBarry Smith   Level: intermediate
7045c6c1daeSBarry Smith 
7055c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
7065c6c1daeSBarry Smith @*/
7075c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
7085c6c1daeSBarry Smith {
7095c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7105c6c1daeSBarry Smith 
7115c6c1daeSBarry Smith   PetscFunctionBegin;
7125c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7135c6c1daeSBarry Smith   hdf5->timestep = timestep;
7145c6c1daeSBarry Smith   PetscFunctionReturn(0);
7155c6c1daeSBarry Smith }
7165c6c1daeSBarry Smith 
7173ef9c667SSatish Balay /*@
7185c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
7195c6c1daeSBarry Smith 
7205c6c1daeSBarry Smith   Not collective
7215c6c1daeSBarry Smith 
7225c6c1daeSBarry Smith   Input Parameter:
7235c6c1daeSBarry Smith . viewer - the PetscViewer
7245c6c1daeSBarry Smith 
7255c6c1daeSBarry Smith   Output Parameter:
7265c6c1daeSBarry Smith . timestep - The timestep number
7275c6c1daeSBarry Smith 
7285c6c1daeSBarry Smith   Level: intermediate
7295c6c1daeSBarry Smith 
7305c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
7315c6c1daeSBarry Smith @*/
7325c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
7335c6c1daeSBarry Smith {
7345c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7355c6c1daeSBarry Smith 
7365c6c1daeSBarry Smith   PetscFunctionBegin;
7375c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7385c6c1daeSBarry Smith   PetscValidPointer(timestep,2);
7395c6c1daeSBarry Smith   *timestep = hdf5->timestep;
7405c6c1daeSBarry Smith   PetscFunctionReturn(0);
7415c6c1daeSBarry Smith }
7425c6c1daeSBarry Smith 
74336bce6e8SMatthew G. Knepley /*@C
74436bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
74536bce6e8SMatthew G. Knepley 
74636bce6e8SMatthew G. Knepley   Not collective
74736bce6e8SMatthew G. Knepley 
74836bce6e8SMatthew G. Knepley   Input Parameter:
74936bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
75036bce6e8SMatthew G. Knepley 
75136bce6e8SMatthew G. Knepley   Output Parameter:
75236bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
75336bce6e8SMatthew G. Knepley 
75436bce6e8SMatthew G. Knepley   Level: advanced
75536bce6e8SMatthew G. Knepley 
75636bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
75736bce6e8SMatthew G. Knepley @*/
75836bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
75936bce6e8SMatthew G. Knepley {
76036bce6e8SMatthew G. Knepley   PetscFunctionBegin;
76136bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
76236bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
76336bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_LLONG;
76436bce6e8SMatthew G. Knepley #else
76536bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_INT;
76636bce6e8SMatthew G. Knepley #endif
76736bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
76836bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
76936bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
77036bce6e8SMatthew G. Knepley   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_DOUBLE;
771de3ba262SVaclav Hapla   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_INT;
77236bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
77336bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
77436bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
7757e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
77636bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
77736bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
77836bce6e8SMatthew G. Knepley }
77936bce6e8SMatthew G. Knepley 
78036bce6e8SMatthew G. Knepley /*@C
78136bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
78236bce6e8SMatthew G. Knepley 
78336bce6e8SMatthew G. Knepley   Not collective
78436bce6e8SMatthew G. Knepley 
78536bce6e8SMatthew G. Knepley   Input Parameter:
78636bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
78736bce6e8SMatthew G. Knepley 
78836bce6e8SMatthew G. Knepley   Output Parameter:
78936bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
79036bce6e8SMatthew G. Knepley 
79136bce6e8SMatthew G. Knepley   Level: advanced
79236bce6e8SMatthew G. Knepley 
79336bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
79436bce6e8SMatthew G. Knepley @*/
79536bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
79636bce6e8SMatthew G. Knepley {
79736bce6e8SMatthew G. Knepley   PetscFunctionBegin;
79836bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
79936bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
80036bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
80136bce6e8SMatthew G. Knepley #else
80236bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
80336bce6e8SMatthew G. Knepley #endif
80436bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
80536bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
80636bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
80736bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
80836bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
80936bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
8107e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
81136bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
81236bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
81336bce6e8SMatthew G. Knepley }
81436bce6e8SMatthew G. Knepley 
815df863907SAlex Fikl /*@C
816b8ee85a1SVaclav Hapla  PetscViewerHDF5WriteAttribute - Write an attribute
81736bce6e8SMatthew G. Knepley 
81836bce6e8SMatthew G. Knepley   Input Parameters:
81936bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer
82057d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
82136bce6e8SMatthew G. Knepley . name   - The attribute name
82236bce6e8SMatthew G. Knepley . datatype - The attribute type
82336bce6e8SMatthew G. Knepley - value    - The attribute value
82436bce6e8SMatthew G. Knepley 
82536bce6e8SMatthew G. Knepley   Level: advanced
82636bce6e8SMatthew G. Knepley 
827577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
82836bce6e8SMatthew G. Knepley @*/
82957d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value)
83036bce6e8SMatthew G. Knepley {
83157d22f7dSVaclav Hapla   char           *parent;
83260568592SMatthew G. Knepley   hid_t          h5, dataspace, obj, attribute, dtype;
833080f144cSVaclav Hapla   PetscBool      has;
83436bce6e8SMatthew G. Knepley   PetscErrorCode ierr;
83536bce6e8SMatthew G. Knepley 
83636bce6e8SMatthew G. Knepley   PetscFunctionBegin;
8375cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
83857d22f7dSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset, 2);
839c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
840b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
84157d22f7dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
842bb286ee1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr);
843b67695ecSVaclav Hapla   ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);
84436bce6e8SMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
8457e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
8467e97c476SMatthew G. Knepley     size_t len;
8477e97c476SMatthew G. Knepley     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
848729ab6d9SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
8497e97c476SMatthew G. Knepley   }
85036bce6e8SMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
851729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
85260568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
853080f144cSVaclav Hapla   if (has) {
854080f144cSVaclav Hapla     PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
855080f144cSVaclav Hapla   } else {
85660568592SMatthew G. Knepley     PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
857080f144cSVaclav Hapla   }
858729ab6d9SBarry Smith   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
859729ab6d9SBarry Smith   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
860729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
86160568592SMatthew G. Knepley   PetscStackCallHDF5(H5Oclose,(obj));
862729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(dataspace));
86357d22f7dSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
86436bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
86536bce6e8SMatthew G. Knepley }
86636bce6e8SMatthew G. Knepley 
867df863907SAlex Fikl /*@C
868577e0e04SVaclav Hapla  PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name
869577e0e04SVaclav Hapla 
870577e0e04SVaclav Hapla   Input Parameters:
871577e0e04SVaclav Hapla + viewer   - The HDF5 viewer
872577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
873577e0e04SVaclav Hapla . name     - The attribute name
874577e0e04SVaclav Hapla . datatype - The attribute type
875577e0e04SVaclav Hapla - value    - The attribute value
876577e0e04SVaclav Hapla 
877577e0e04SVaclav Hapla   Notes:
878577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
879577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
880577e0e04SVaclav Hapla 
881577e0e04SVaclav Hapla   Level: advanced
882577e0e04SVaclav Hapla 
883577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
884577e0e04SVaclav Hapla @*/
885577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
886577e0e04SVaclav Hapla {
887577e0e04SVaclav Hapla   PetscErrorCode ierr;
888577e0e04SVaclav Hapla 
889577e0e04SVaclav Hapla   PetscFunctionBegin;
890577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
891577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
892577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
893b7e81f0fSVaclav Hapla   PetscValidPointer(value,5);
894577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
895577e0e04SVaclav Hapla   ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
896577e0e04SVaclav Hapla   PetscFunctionReturn(0);
897577e0e04SVaclav Hapla }
898577e0e04SVaclav Hapla 
899577e0e04SVaclav Hapla /*@C
900b8ee85a1SVaclav Hapla  PetscViewerHDF5ReadAttribute - Read an attribute
901ad0c61feSMatthew G. Knepley 
902ad0c61feSMatthew G. Knepley   Input Parameters:
903ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer
90457d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
905ad0c61feSMatthew G. Knepley . name   - The attribute name
906ad0c61feSMatthew G. Knepley - datatype - The attribute type
907ad0c61feSMatthew G. Knepley 
908ad0c61feSMatthew G. Knepley   Output Parameter:
909ad0c61feSMatthew G. Knepley . value    - The attribute value
910ad0c61feSMatthew G. Knepley 
911ad0c61feSMatthew G. Knepley   Level: advanced
912ad0c61feSMatthew G. Knepley 
913577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
914ad0c61feSMatthew G. Knepley @*/
91557d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value)
916ad0c61feSMatthew G. Knepley {
91757d22f7dSVaclav Hapla   char           *parent;
918f0b58479SMatthew G. Knepley   hid_t          h5, obj, attribute, atype, dtype;
919969748fdSVaclav Hapla   PetscBool      has;
920ad0c61feSMatthew G. Knepley   PetscErrorCode ierr;
921ad0c61feSMatthew G. Knepley 
922ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
9235cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
92457d22f7dSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset, 2);
925c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
926b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
92757d22f7dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
928969748fdSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr);
929969748fdSVaclav Hapla   if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);}
930969748fdSVaclav Hapla   if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name);
931ad0c61feSMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
932ad0c61feSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
93360568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
93460568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
935f0b58479SMatthew G. Knepley   if (datatype == PETSC_STRING) {
936f0b58479SMatthew G. Knepley     size_t len;
937e90c5075SMarek Pecha     PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
938f0b58479SMatthew G. Knepley     PetscStackCallHDF5Return(len,H5Tget_size,(atype));
939f0b58479SMatthew G. Knepley     PetscStackCallHDF5(H5Tclose,(atype));
940f0b58479SMatthew G. Knepley     ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr);
941f0b58479SMatthew G. Knepley   }
94270efba33SBarry Smith   PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
943729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
944e90c5075SMarek Pecha   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
945e90c5075SMarek Pecha   PetscStackCallHDF5(H5Oclose,(obj));
94657d22f7dSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
947ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
948ad0c61feSMatthew G. Knepley }
949ad0c61feSMatthew G. Knepley 
950577e0e04SVaclav Hapla /*@C
951577e0e04SVaclav Hapla  PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name
952577e0e04SVaclav Hapla 
953577e0e04SVaclav Hapla   Input Parameters:
954577e0e04SVaclav Hapla + viewer   - The HDF5 viewer
955577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
956577e0e04SVaclav Hapla . name     - The attribute name
957577e0e04SVaclav Hapla - datatype - The attribute type
958577e0e04SVaclav Hapla 
959577e0e04SVaclav Hapla   Output Parameter:
960577e0e04SVaclav Hapla . value    - The attribute value
961577e0e04SVaclav Hapla 
962577e0e04SVaclav Hapla   Notes:
963577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
964577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
965577e0e04SVaclav Hapla 
966577e0e04SVaclav Hapla   Level: advanced
967577e0e04SVaclav Hapla 
968577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
969577e0e04SVaclav Hapla @*/
970577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value)
971577e0e04SVaclav Hapla {
972577e0e04SVaclav Hapla   PetscErrorCode ierr;
973577e0e04SVaclav Hapla 
974577e0e04SVaclav Hapla   PetscFunctionBegin;
975577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
976577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
977577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
978b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
979577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
980577e0e04SVaclav Hapla   ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
981577e0e04SVaclav Hapla   PetscFunctionReturn(0);
982577e0e04SVaclav Hapla }
983577e0e04SVaclav Hapla 
984a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
985a75c3fd4SVaclav Hapla {
986a75c3fd4SVaclav Hapla   htri_t exists;
987a75c3fd4SVaclav Hapla   hid_t group;
988a75c3fd4SVaclav Hapla 
989a75c3fd4SVaclav Hapla   PetscFunctionBegin;
990a75c3fd4SVaclav Hapla   PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT));
991a75c3fd4SVaclav Hapla   if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT));
992a75c3fd4SVaclav Hapla   if (!exists && createGroup) {
993a75c3fd4SVaclav Hapla     PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
994a75c3fd4SVaclav Hapla     PetscStackCallHDF5(H5Gclose,(group));
995a75c3fd4SVaclav Hapla     exists = PETSC_TRUE;
996a75c3fd4SVaclav Hapla   }
997a75c3fd4SVaclav Hapla   *exists_ = (PetscBool) exists;
998a75c3fd4SVaclav Hapla   PetscFunctionReturn(0);
999a75c3fd4SVaclav Hapla }
1000a75c3fd4SVaclav Hapla 
1001bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
10025cdeb108SMatthew G. Knepley {
100390e03537SVaclav Hapla   const char     rootGroupName[] = "/";
10045cdeb108SMatthew G. Knepley   hid_t          h5;
1005e5bf9ebcSVaclav Hapla   PetscBool      exists=PETSC_FALSE;
100685688ae2SVaclav Hapla   PetscInt       i,n;
100785688ae2SVaclav Hapla   char           **hierarchy;
100885688ae2SVaclav Hapla   char           buf[PETSC_MAX_PATH_LEN]="";
10095cdeb108SMatthew G. Knepley   PetscErrorCode ierr;
10105cdeb108SMatthew G. Knepley 
10115cdeb108SMatthew G. Knepley   PetscFunctionBegin;
10125cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
101390e03537SVaclav Hapla   if (name) PetscValidCharPointer(name, 2);
101490e03537SVaclav Hapla   else name = rootGroupName;
1015ccd66a83SVaclav Hapla   if (has) {
101656cc0592SVaclav Hapla     PetscValidIntPointer(has, 3);
10175cdeb108SMatthew G. Knepley     *has = PETSC_FALSE;
1018ccd66a83SVaclav Hapla   }
1019ccd66a83SVaclav Hapla   if (otype) {
1020ccd66a83SVaclav Hapla     PetscValidIntPointer(otype, 4);
102156cc0592SVaclav Hapla     *otype = H5O_TYPE_UNKNOWN;
1022ccd66a83SVaclav Hapla   }
1023c0cd0301SJed Brown   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
102485688ae2SVaclav Hapla 
102585688ae2SVaclav Hapla   /*
102685688ae2SVaclav Hapla      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
102785688ae2SVaclav Hapla      Hence, each of them needs to be tested separately:
102885688ae2SVaclav Hapla      1) whether it's a valid link
102985688ae2SVaclav Hapla      2) whether this link resolves to an object
103085688ae2SVaclav Hapla      See H5Oexists_by_name() documentation.
103185688ae2SVaclav Hapla   */
103285688ae2SVaclav Hapla   ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr);
103385688ae2SVaclav Hapla   if (!n) {
103485688ae2SVaclav Hapla     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1035ccd66a83SVaclav Hapla     if (has)   *has   = PETSC_TRUE;
1036ccd66a83SVaclav Hapla     if (otype) *otype = H5O_TYPE_GROUP;
103785688ae2SVaclav Hapla     ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
103885688ae2SVaclav Hapla     PetscFunctionReturn(0);
103985688ae2SVaclav Hapla   }
104085688ae2SVaclav Hapla   for (i=0; i<n; i++) {
104185688ae2SVaclav Hapla     ierr = PetscStrcat(buf,"/");CHKERRQ(ierr);
104285688ae2SVaclav Hapla     ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr);
1043a75c3fd4SVaclav Hapla     ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr);
1044a75c3fd4SVaclav Hapla     if (!exists) break;
104585688ae2SVaclav Hapla   }
104685688ae2SVaclav Hapla   ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
104785688ae2SVaclav Hapla 
104885688ae2SVaclav Hapla   /* If the object exists, get its type */
1049ccd66a83SVaclav Hapla   if (exists && otype) {
10505cdeb108SMatthew G. Knepley     H5O_info_t info;
10515cdeb108SMatthew G. Knepley 
105276276748SVaclav Hapla     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
105304633f7fSVaclav Hapla     PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT));
105456cc0592SVaclav Hapla     *otype = info.type;
10555cdeb108SMatthew G. Knepley   }
1056ccd66a83SVaclav Hapla   if (has) *has = exists;
10575cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
10585cdeb108SMatthew G. Knepley }
10595cdeb108SMatthew G. Knepley 
1060ecce1506SVaclav Hapla /*@
10610a9f38efSVaclav Hapla  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
10620a9f38efSVaclav Hapla 
10630a9f38efSVaclav Hapla   Input Parameters:
10640a9f38efSVaclav Hapla . viewer - The HDF5 viewer
10650a9f38efSVaclav Hapla 
10660a9f38efSVaclav Hapla   Output Parameter:
10670a9f38efSVaclav Hapla . has    - Flag for group existence
10680a9f38efSVaclav Hapla 
10690a9f38efSVaclav Hapla   Notes:
10700a9f38efSVaclav Hapla   If the path exists but is not a group, this returns PETSC_FALSE as well.
10710a9f38efSVaclav Hapla 
10720a9f38efSVaclav Hapla   Level: advanced
10730a9f38efSVaclav Hapla 
10740a9f38efSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup()
10750a9f38efSVaclav Hapla @*/
10760a9f38efSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has)
10770a9f38efSVaclav Hapla {
10780a9f38efSVaclav Hapla   H5O_type_t type;
10790a9f38efSVaclav Hapla   const char *name;
10800a9f38efSVaclav Hapla   PetscErrorCode ierr;
10810a9f38efSVaclav Hapla 
10820a9f38efSVaclav Hapla   PetscFunctionBegin;
10830a9f38efSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
10840a9f38efSVaclav Hapla   PetscValidIntPointer(has,2);
10850a9f38efSVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr);
10860a9f38efSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr);
10870a9f38efSVaclav Hapla   *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE;
10880a9f38efSVaclav Hapla   PetscFunctionReturn(0);
10890a9f38efSVaclav Hapla }
10900a9f38efSVaclav Hapla 
10910a9f38efSVaclav Hapla /*@
1092e3f143f7SVaclav Hapla  PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1093ecce1506SVaclav Hapla 
1094ecce1506SVaclav Hapla   Input Parameters:
1095ecce1506SVaclav Hapla + viewer - The HDF5 viewer
1096ecce1506SVaclav Hapla - obj    - The named object
1097ecce1506SVaclav Hapla 
1098ecce1506SVaclav Hapla   Output Parameter:
1099ecce1506SVaclav Hapla . has    - Flag for dataset existence; PETSC_FALSE for unnamed object
1100ecce1506SVaclav Hapla 
1101e3f143f7SVaclav Hapla   Notes:
1102e3f143f7SVaclav Hapla   If the path exists but is not a dataset, this returns PETSC_FALSE as well.
1103e3f143f7SVaclav Hapla 
1104ecce1506SVaclav Hapla   Level: advanced
1105ecce1506SVaclav Hapla 
1106e3f143f7SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1107ecce1506SVaclav Hapla @*/
1108ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1109ecce1506SVaclav Hapla {
111056cc0592SVaclav Hapla   H5O_type_t type;
11116c132bc1SVaclav Hapla   char *path;
1112ecce1506SVaclav Hapla   PetscErrorCode ierr;
1113ecce1506SVaclav Hapla 
1114ecce1506SVaclav Hapla   PetscFunctionBegin;
1115c57a1a9aSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1116c57a1a9aSVaclav Hapla   PetscValidHeader(obj,2);
1117c57a1a9aSVaclav Hapla   PetscValidIntPointer(has,3);
1118ecce1506SVaclav Hapla   *has = PETSC_FALSE;
1119e3f143f7SVaclav Hapla   if (!obj->name) PetscFunctionReturn(0);
11206c132bc1SVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr);
11216c132bc1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr);
112256cc0592SVaclav Hapla   *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE;
11236c132bc1SVaclav Hapla   ierr = PetscFree(path);CHKERRQ(ierr);
1124ecce1506SVaclav Hapla   PetscFunctionReturn(0);
1125ecce1506SVaclav Hapla }
1126ecce1506SVaclav Hapla 
1127df863907SAlex Fikl /*@C
1128ece83b6cSVaclav Hapla  PetscViewerHDF5HasAttribute - Check whether an attribute exists
1129e7bdbf8eSMatthew G. Knepley 
1130e7bdbf8eSMatthew G. Knepley   Input Parameters:
1131e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer
113210e69e8fSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
1133e7bdbf8eSMatthew G. Knepley - name   - The attribute name
1134e7bdbf8eSMatthew G. Knepley 
1135e7bdbf8eSMatthew G. Knepley   Output Parameter:
1136e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
1137e7bdbf8eSMatthew G. Knepley 
1138e7bdbf8eSMatthew G. Knepley   Level: advanced
1139e7bdbf8eSMatthew G. Knepley 
1140577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1141e7bdbf8eSMatthew G. Knepley @*/
114210e69e8fSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has)
1143e7bdbf8eSMatthew G. Knepley {
114410e69e8fSVaclav Hapla   char           *parent;
1145e7bdbf8eSMatthew G. Knepley   PetscErrorCode ierr;
1146e7bdbf8eSMatthew G. Knepley 
1147e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
11485cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
114910e69e8fSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset,2);
1150c57a1a9aSVaclav Hapla   PetscValidCharPointer(name,3);
1151c57a1a9aSVaclav Hapla   PetscValidIntPointer(has,4);
115210e69e8fSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
1153bb286ee1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr);
115410e69e8fSVaclav Hapla   if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);}
115510e69e8fSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
115606db490cSVaclav Hapla   PetscFunctionReturn(0);
115706db490cSVaclav Hapla }
115806db490cSVaclav Hapla 
1159577e0e04SVaclav Hapla /*@C
1160577e0e04SVaclav Hapla  PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name
1161577e0e04SVaclav Hapla 
1162577e0e04SVaclav Hapla   Input Parameters:
1163577e0e04SVaclav Hapla + viewer - The HDF5 viewer
1164577e0e04SVaclav Hapla . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1165577e0e04SVaclav Hapla - name   - The attribute name
1166577e0e04SVaclav Hapla 
1167577e0e04SVaclav Hapla   Output Parameter:
1168577e0e04SVaclav Hapla . has    - Flag for attribute existence
1169577e0e04SVaclav Hapla 
1170577e0e04SVaclav Hapla   Notes:
1171577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1172577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
1173577e0e04SVaclav Hapla 
1174577e0e04SVaclav Hapla   Level: advanced
1175577e0e04SVaclav Hapla 
1176577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1177577e0e04SVaclav Hapla @*/
1178577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1179577e0e04SVaclav Hapla {
1180577e0e04SVaclav Hapla   PetscErrorCode ierr;
1181577e0e04SVaclav Hapla 
1182577e0e04SVaclav Hapla   PetscFunctionBegin;
1183577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1184577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
1185577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
1186577e0e04SVaclav Hapla   PetscValidIntPointer(has,4);
1187577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
1188577e0e04SVaclav Hapla   ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr);
1189577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1190577e0e04SVaclav Hapla }
1191577e0e04SVaclav Hapla 
119206db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
119306db490cSVaclav Hapla {
119406db490cSVaclav Hapla   hid_t          h5;
119506db490cSVaclav Hapla   htri_t         hhas;
119606db490cSVaclav Hapla   PetscErrorCode ierr;
119706db490cSVaclav Hapla 
119806db490cSVaclav Hapla   PetscFunctionBegin;
119906db490cSVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
12002f936e54SVaclav Hapla   PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT));
12015cdeb108SMatthew G. Knepley   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1202e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
1203e7bdbf8eSMatthew G. Knepley }
1204e7bdbf8eSMatthew G. Knepley 
1205eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadInitialize_Private(PetscViewer viewer, const char name[], HDF5ReadCtx *ctx)
1206a5e1feadSVaclav Hapla {
1207fbd37863SVaclav Hapla   HDF5ReadCtx    h=NULL;
1208a5e1feadSVaclav Hapla   PetscErrorCode ierr;
1209a5e1feadSVaclav Hapla 
1210a5e1feadSVaclav Hapla   PetscFunctionBegin;
121169a06e7bSVaclav Hapla   ierr = PetscNew(&h);CHKERRQ(ierr);
121269a06e7bSVaclav Hapla   ierr = PetscViewerHDF5OpenGroup(viewer, &h->file, &h->group);CHKERRQ(ierr);
121369a06e7bSVaclav Hapla   PetscStackCallHDF5Return(h->dataset,H5Dopen2,(h->group, name, H5P_DEFAULT));
121469a06e7bSVaclav Hapla   PetscStackCallHDF5Return(h->dataspace,H5Dget_space,(h->dataset));
12159568d583SVaclav Hapla   ierr = PetscViewerHDF5GetTimestep(viewer, &h->timestep);CHKERRQ(ierr);
1216adb363a2SVaclav Hapla   ierr = PetscViewerHDF5HasAttribute(viewer,name,"complex",&h->complexVal);CHKERRQ(ierr);
1217aa54fc5aSVaclav Hapla   if (h->complexVal) {ierr = PetscViewerHDF5ReadAttribute(viewer,name,"complex",PETSC_BOOL,&h->complexVal);CHKERRQ(ierr);}
121809dabeb0SVaclav Hapla   /* MATLAB stores column vectors horizontally */
1219adb363a2SVaclav Hapla   ierr = PetscViewerHDF5HasAttribute(viewer,name,"MATLAB_class",&h->horizontal);CHKERRQ(ierr);
1220b8ef406cSVaclav Hapla   /* Create property list for collective dataset read */
1221b8ef406cSVaclav Hapla   PetscStackCallHDF5Return(h->plist,H5Pcreate,(H5P_DATASET_XFER));
1222b8ef406cSVaclav Hapla #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
1223b8ef406cSVaclav Hapla   PetscStackCallHDF5(H5Pset_dxpl_mpio,(h->plist, H5FD_MPIO_COLLECTIVE));
1224b8ef406cSVaclav Hapla #endif
1225b8ef406cSVaclav Hapla   /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */
122669a06e7bSVaclav Hapla   *ctx = h;
122769a06e7bSVaclav Hapla   PetscFunctionReturn(0);
122869a06e7bSVaclav Hapla }
122969a06e7bSVaclav Hapla 
1230eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadFinalize_Private(PetscViewer viewer, HDF5ReadCtx *ctx)
123169a06e7bSVaclav Hapla {
123269a06e7bSVaclav Hapla   HDF5ReadCtx    h;
123369a06e7bSVaclav Hapla   PetscErrorCode ierr;
123469a06e7bSVaclav Hapla 
123569a06e7bSVaclav Hapla   PetscFunctionBegin;
123669a06e7bSVaclav Hapla   h = *ctx;
1237b8ef406cSVaclav Hapla   PetscStackCallHDF5(H5Pclose,(h->plist));
123890e03537SVaclav Hapla   PetscStackCallHDF5(H5Gclose,(h->group));
123969a06e7bSVaclav Hapla   PetscStackCallHDF5(H5Sclose,(h->dataspace));
124069a06e7bSVaclav Hapla   PetscStackCallHDF5(H5Dclose,(h->dataset));
124169a06e7bSVaclav Hapla   ierr = PetscFree(*ctx);CHKERRQ(ierr);
124269a06e7bSVaclav Hapla   PetscFunctionReturn(0);
124369a06e7bSVaclav Hapla }
124469a06e7bSVaclav Hapla 
1245eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadSizes_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout *map_)
124669a06e7bSVaclav Hapla {
124769a06e7bSVaclav Hapla   int            rdim, dim;
124869a06e7bSVaclav Hapla   hsize_t        dims[4];
124909dabeb0SVaclav Hapla   PetscInt       bsInd, lenInd, bs, len, N;
12508374c777SVaclav Hapla   PetscLayout    map;
12518374c777SVaclav Hapla   PetscErrorCode ierr;
125269a06e7bSVaclav Hapla 
125369a06e7bSVaclav Hapla   PetscFunctionBegin;
12548374c777SVaclav Hapla   if (!(*map_)) {
12558374c777SVaclav Hapla     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)viewer),map_);CHKERRQ(ierr);
12568374c777SVaclav Hapla   }
12578374c777SVaclav Hapla   map = *map_;
12589787f08bSVaclav Hapla   /* calculate expected number of dimensions */
1259a5e1feadSVaclav Hapla   dim = 0;
12609568d583SVaclav Hapla   if (ctx->timestep >= 0) ++dim;
1261a5e1feadSVaclav Hapla   ++dim; /* length in blocks */
12629568d583SVaclav Hapla   if (ctx->complexVal) ++dim;
12639787f08bSVaclav Hapla   /* get actual number of dimensions in dataset */
126469a06e7bSVaclav Hapla   PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(ctx->dataspace, dims, NULL));
12659787f08bSVaclav Hapla   /* calculate expected dimension indices */
12669787f08bSVaclav Hapla   lenInd = 0;
12679568d583SVaclav Hapla   if (ctx->timestep >= 0) ++lenInd;
12689787f08bSVaclav Hapla   bsInd = lenInd + 1;
12699568d583SVaclav Hapla   ctx->dim2 = PETSC_FALSE;
12709787f08bSVaclav Hapla   if (rdim == dim) {
127145e21b7fSVaclav Hapla     bs = 1; /* support vectors stored as 1D array */
12729787f08bSVaclav Hapla   } else if (rdim == dim+1) {
127345e21b7fSVaclav Hapla     bs = (PetscInt) dims[bsInd];
12749568d583SVaclav Hapla     if (bs == 1) ctx->dim2 = PETSC_TRUE; /* vector with blocksize of 1, still stored as 2D array */
1275a5e1feadSVaclav Hapla   } else {
1276b54f1188SVaclav Hapla     SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Number of dimensions %d not %d as expected", rdim, dim);
1277a5e1feadSVaclav Hapla   }
127809dabeb0SVaclav Hapla   len = dims[lenInd];
127909dabeb0SVaclav Hapla   if (ctx->horizontal) {
1280b54f1188SVaclav Hapla     if (len != 1) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Cannot have horizontal array with number of rows > 1. In case of MATLAB MAT-file, vectors must be saved as column vectors.");
128109dabeb0SVaclav Hapla     len = bs;
128209dabeb0SVaclav Hapla     bs = 1;
128309dabeb0SVaclav Hapla   }
128409dabeb0SVaclav Hapla   N = (PetscInt) len*bs;
12858374c777SVaclav Hapla 
12868374c777SVaclav Hapla   /* Set Vec sizes,blocksize,and type if not already set */
12878374c777SVaclav Hapla   if (map->bs < 0) {
128845e21b7fSVaclav Hapla     ierr = PetscLayoutSetBlockSize(map, bs);CHKERRQ(ierr);
128945e21b7fSVaclav 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);
12908374c777SVaclav Hapla   if (map->N < 0) {
129145e21b7fSVaclav Hapla     ierr = PetscLayoutSetSize(map, N);CHKERRQ(ierr);
1292b54f1188SVaclav Hapla   } else if (map->N != N) SETERRQ2(PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED, "Global size of array in file is %D, not %D as expected",N,map->N);
129369a06e7bSVaclav Hapla   PetscFunctionReturn(0);
129469a06e7bSVaclav Hapla }
129569a06e7bSVaclav Hapla 
1296eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadSelectHyperslab_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout map, hid_t *memspace)
129716f6402dSVaclav Hapla {
129816f6402dSVaclav Hapla   hsize_t        count[4], offset[4];
129916f6402dSVaclav Hapla   int            dim;
130016f6402dSVaclav Hapla   PetscInt       bs, n, low;
130116f6402dSVaclav Hapla   PetscErrorCode ierr;
130216f6402dSVaclav Hapla 
130316f6402dSVaclav Hapla   PetscFunctionBegin;
130416f6402dSVaclav Hapla   /* Compute local size and ownership range */
130516f6402dSVaclav Hapla   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
130616f6402dSVaclav Hapla   ierr = PetscLayoutGetBlockSize(map, &bs);CHKERRQ(ierr);
130716f6402dSVaclav Hapla   ierr = PetscLayoutGetLocalSize(map, &n);CHKERRQ(ierr);
130816f6402dSVaclav Hapla   ierr = PetscLayoutGetRange(map, &low, NULL);CHKERRQ(ierr);
130916f6402dSVaclav Hapla 
131016f6402dSVaclav Hapla   /* Each process defines a dataset and reads it from the hyperslab in the file */
131116f6402dSVaclav Hapla   dim  = 0;
131216f6402dSVaclav Hapla   if (ctx->timestep >= 0) {
131316f6402dSVaclav Hapla     count[dim]  = 1;
131416f6402dSVaclav Hapla     offset[dim] = ctx->timestep;
131516f6402dSVaclav Hapla     ++dim;
131616f6402dSVaclav Hapla   }
131709dabeb0SVaclav Hapla   if (ctx->horizontal) {
131809dabeb0SVaclav Hapla     count[dim]  = 1;
131909dabeb0SVaclav Hapla     offset[dim] = 0;
132009dabeb0SVaclav Hapla     ++dim;
132109dabeb0SVaclav Hapla   }
132216f6402dSVaclav Hapla   {
132316f6402dSVaclav Hapla     ierr = PetscHDF5IntCast(n/bs, &count[dim]);CHKERRQ(ierr);
132416f6402dSVaclav Hapla     ierr = PetscHDF5IntCast(low/bs, &offset[dim]);CHKERRQ(ierr);
132516f6402dSVaclav Hapla     ++dim;
132616f6402dSVaclav Hapla   }
132716f6402dSVaclav Hapla   if (bs > 1 || ctx->dim2) {
132809dabeb0SVaclav Hapla     if (PetscUnlikely(ctx->horizontal)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "cannot have horizontal array with blocksize > 1");
132916f6402dSVaclav Hapla     count[dim]  = bs;
133016f6402dSVaclav Hapla     offset[dim] = 0;
133116f6402dSVaclav Hapla     ++dim;
133216f6402dSVaclav Hapla   }
133316f6402dSVaclav Hapla   if (ctx->complexVal) {
133416f6402dSVaclav Hapla     count[dim]  = 2;
133516f6402dSVaclav Hapla     offset[dim] = 0;
133616f6402dSVaclav Hapla     ++dim;
133716f6402dSVaclav Hapla   }
133816f6402dSVaclav Hapla   PetscStackCallHDF5Return(*memspace,H5Screate_simple,(dim, count, NULL));
133916f6402dSVaclav Hapla   PetscStackCallHDF5(H5Sselect_hyperslab,(ctx->dataspace, H5S_SELECT_SET, offset, NULL, count, NULL));
134016f6402dSVaclav Hapla   PetscFunctionReturn(0);
134116f6402dSVaclav Hapla }
134216f6402dSVaclav Hapla 
1343eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadArray_Private(PetscViewer viewer, HDF5ReadCtx h, hid_t datatype, hid_t memspace, void *arr)
1344ef2d82ceSVaclav Hapla {
1345ef2d82ceSVaclav Hapla   PetscFunctionBegin;
1346ef2d82ceSVaclav Hapla   PetscStackCallHDF5(H5Dread,(h->dataset, datatype, memspace, h->dataspace, h->plist, arr));
1347ef2d82ceSVaclav Hapla   PetscFunctionReturn(0);
1348ef2d82ceSVaclav Hapla }
1349ef2d82ceSVaclav Hapla 
1350dad982a8SVaclav Hapla PetscErrorCode PetscViewerHDF5Load(PetscViewer viewer, const char *name, PetscLayout map, hid_t datatype, void **newarr)
1351a25c73c6SVaclav Hapla {
1352fbd37863SVaclav Hapla   HDF5ReadCtx     h=NULL;
1353fbd37863SVaclav Hapla   hid_t           memspace=0;
1354a25c73c6SVaclav Hapla   size_t          unitsize;
1355a25c73c6SVaclav Hapla   void            *arr;
1356a25c73c6SVaclav Hapla   PetscErrorCode  ierr;
1357a25c73c6SVaclav Hapla 
1358a25c73c6SVaclav Hapla   PetscFunctionBegin;
1359eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr);
13605a89fdf4SVaclav Hapla #if defined(PETSC_USE_COMPLEX)
13615a89fdf4SVaclav Hapla   if (!h->complexVal) {
1362c76551afSVaclav Hapla     H5T_class_t clazz = H5Tget_class(datatype);
1363c76551afSVaclav 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.");
13645a89fdf4SVaclav Hapla   }
13655a89fdf4SVaclav Hapla #else
13665a89fdf4SVaclav 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.");
13675a89fdf4SVaclav Hapla #endif
1368e9e90110SVaclav Hapla 
1369e9e90110SVaclav Hapla   ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr);
1370e9e90110SVaclav Hapla   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
1371e9e90110SVaclav Hapla   ierr = PetscViewerHDF5ReadSelectHyperslab_Private(viewer, h, map, &memspace);CHKERRQ(ierr);
1372e9e90110SVaclav Hapla 
13734fc17bcdSVaclav Hapla   unitsize = H5Tget_size(datatype);
13744fc17bcdSVaclav Hapla   if (h->complexVal) unitsize *= 2;
1375dff35581SVaclav Hapla   if (unitsize <= 0 || unitsize > PetscMax(sizeof(PetscInt),sizeof(PetscScalar))) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Sanity check failed: HDF5 function H5Tget_size(datatype) returned suspicious value %D",unitsize);
13764fc17bcdSVaclav Hapla   ierr = PetscMalloc(map->n*unitsize, &arr);CHKERRQ(ierr);
13774fc17bcdSVaclav Hapla 
1378eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadArray_Private(viewer, h, datatype, memspace, arr);CHKERRQ(ierr);
1379a25c73c6SVaclav Hapla   PetscStackCallHDF5(H5Sclose,(memspace));
1380eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr);
1381a25c73c6SVaclav Hapla   *newarr = arr;
1382a25c73c6SVaclav Hapla   PetscFunctionReturn(0);
1383a25c73c6SVaclav Hapla }
1384a25c73c6SVaclav Hapla 
1385c1aaad9cSVaclav Hapla /*@C
1386c1aaad9cSVaclav Hapla  PetscViewerHDF5ReadSizes - Read block size and global size of a vector (Vec or IS) stored in an HDF5 file.
1387c1aaad9cSVaclav Hapla 
1388c1aaad9cSVaclav Hapla   Input Parameters:
1389c1aaad9cSVaclav Hapla + viewer - The HDF5 viewer
1390c1aaad9cSVaclav Hapla - name   - The vector name
1391c1aaad9cSVaclav Hapla 
1392c1aaad9cSVaclav Hapla   Output Parameter:
1393c1aaad9cSVaclav Hapla + bs     - block size
1394c1aaad9cSVaclav Hapla - N      - global size
1395c1aaad9cSVaclav Hapla 
1396c1aaad9cSVaclav Hapla   Note:
1397c1aaad9cSVaclav Hapla   A vector is stored as an HDF5 dataspace with 1-4 dimensions in this order:
1398c1aaad9cSVaclav Hapla   1) # timesteps (optional), 2) # blocks, 3) # elements per block (optional), 4) real and imaginary part (only for complex).
1399c1aaad9cSVaclav Hapla 
1400c1aaad9cSVaclav Hapla   A vectors can be stored as a 2D dataspace even if its blocksize is 1; see PetscViewerHDF5SetBaseDimension2().
1401c1aaad9cSVaclav Hapla 
1402c1aaad9cSVaclav Hapla   Level: advanced
1403c1aaad9cSVaclav Hapla 
1404c1aaad9cSVaclav Hapla .seealso: PetscViewerHDF5Open(), VecLoad(), ISLoad(), VecGetSize(), ISGetSize(), PetscViewerHDF5SetBaseDimension2()
1405c1aaad9cSVaclav Hapla @*/
140669a06e7bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadSizes(PetscViewer viewer, const char name[], PetscInt *bs, PetscInt *N)
140769a06e7bSVaclav Hapla {
1408fbd37863SVaclav Hapla   HDF5ReadCtx    h=NULL;
14098374c777SVaclav Hapla   PetscLayout    map=NULL;
141069a06e7bSVaclav Hapla   PetscErrorCode ierr;
141169a06e7bSVaclav Hapla 
141269a06e7bSVaclav Hapla   PetscFunctionBegin;
1413c1aaad9cSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1414eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr);
1415eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr);
1416eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr);
14178374c777SVaclav Hapla   if (bs) *bs = map->bs;
14188374c777SVaclav Hapla   if (N) *N = map->N;
14198374c777SVaclav Hapla   ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr);
1420a5e1feadSVaclav Hapla   PetscFunctionReturn(0);
1421a5e1feadSVaclav Hapla }
1422a5e1feadSVaclav Hapla 
1423a75e6a4aSMatthew G. Knepley /*
1424a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1425a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
1426a75e6a4aSMatthew G. Knepley */
1427d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1428a75e6a4aSMatthew G. Knepley 
1429a75e6a4aSMatthew G. Knepley /*@C
1430a75e6a4aSMatthew G. Knepley   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
1431a75e6a4aSMatthew G. Knepley 
1432a75e6a4aSMatthew G. Knepley   Collective on MPI_Comm
1433a75e6a4aSMatthew G. Knepley 
1434a75e6a4aSMatthew G. Knepley   Input Parameter:
1435a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer
1436a75e6a4aSMatthew G. Knepley 
1437a75e6a4aSMatthew G. Knepley   Level: intermediate
1438a75e6a4aSMatthew G. Knepley 
1439a75e6a4aSMatthew G. Knepley   Options Database Keys:
1440a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name>
1441a75e6a4aSMatthew G. Knepley 
1442a75e6a4aSMatthew G. Knepley   Environmental variables:
1443a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME
1444a75e6a4aSMatthew G. Knepley 
1445a75e6a4aSMatthew G. Knepley   Notes:
1446a75e6a4aSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
1447a75e6a4aSMatthew G. Knepley   an error code.  The HDF5 PetscViewer is usually used in the form
1448a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1449a75e6a4aSMatthew G. Knepley 
1450a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1451a75e6a4aSMatthew G. Knepley @*/
1452a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1453a75e6a4aSMatthew G. Knepley {
1454a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
1455a75e6a4aSMatthew G. Knepley   PetscBool      flg;
1456a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
1457a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
1458a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
1459a75e6a4aSMatthew G. Knepley 
1460a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
1461a75e6a4aSMatthew 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);}
1462a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
146312801b39SBarry Smith     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
1464a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1465a75e6a4aSMatthew G. Knepley   }
146647435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1467a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1468a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
1469a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
1470a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1471a75e6a4aSMatthew G. Knepley     if (!flg) {
1472a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname,"output.h5");
1473a75e6a4aSMatthew G. Knepley       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1474a75e6a4aSMatthew G. Knepley     }
1475a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
1476a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1477a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
1478a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
147947435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
1480a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1481a75e6a4aSMatthew G. Knepley   }
1482a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
1483a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1484a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
1485a75e6a4aSMatthew G. Knepley }
1486