xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 6b9fdc0ace1199a7b9ffdadd61bc62ae7357afe7)
1*6b9fdc0aSVaclav Hapla #include <petsc/private/viewerhdf5impl.h> /*I "petscviewerhdf5.h" I*/
25c6c1daeSBarry Smith 
3bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool*, H5O_type_t*);
406db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool*);
506db490cSVaclav Hapla 
64302210dSVaclav Hapla static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char path[], char *abspath[])
76c132bc1SVaclav Hapla {
84302210dSVaclav Hapla   PetscBool      relative = PETSC_FALSE;
96c132bc1SVaclav Hapla   const char     *group;
106c132bc1SVaclav Hapla   char           buf[PETSC_MAX_PATH_LEN] = "";
116c132bc1SVaclav Hapla   PetscErrorCode ierr;
126c132bc1SVaclav Hapla 
136c132bc1SVaclav Hapla   PetscFunctionBegin;
146c132bc1SVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr);
154302210dSVaclav Hapla   ierr = PetscViewerHDF5PathIsRelative(path, PETSC_TRUE, &relative);CHKERRQ(ierr);
164302210dSVaclav Hapla   if (relative) {
174302210dSVaclav Hapla     ierr = PetscStrcpy(buf, group);CHKERRQ(ierr);
186c132bc1SVaclav Hapla     ierr = PetscStrcat(buf, "/");CHKERRQ(ierr);
194302210dSVaclav Hapla     ierr = PetscStrcat(buf, path);CHKERRQ(ierr);
204302210dSVaclav Hapla     ierr = PetscStrallocpy(buf, abspath);CHKERRQ(ierr);
214302210dSVaclav Hapla   } else {
224302210dSVaclav Hapla     ierr = PetscStrallocpy(path, abspath);CHKERRQ(ierr);
234302210dSVaclav Hapla   }
246c132bc1SVaclav Hapla   PetscFunctionReturn(0);
256c132bc1SVaclav Hapla }
266c132bc1SVaclav Hapla 
27577e0e04SVaclav Hapla static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj)
28577e0e04SVaclav Hapla {
29577e0e04SVaclav Hapla   PetscBool      has;
30577e0e04SVaclav Hapla   PetscErrorCode ierr;
31577e0e04SVaclav Hapla 
32577e0e04SVaclav Hapla   PetscFunctionBegin;
33577e0e04SVaclav Hapla   ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr);
34577e0e04SVaclav Hapla   if (!has) {
3589e0ef10SVaclav Hapla     const char *group;
36577e0e04SVaclav Hapla     ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr);
3789e0ef10SVaclav Hapla     SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", obj->name, group ? group : "/");
38577e0e04SVaclav Hapla   }
39577e0e04SVaclav Hapla   PetscFunctionReturn(0);
40577e0e04SVaclav Hapla }
41577e0e04SVaclav Hapla 
424416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v)
4382ea9b62SBarry Smith {
4482ea9b62SBarry Smith   PetscErrorCode   ierr;
45ee8b9732SVaclav Hapla   PetscBool        flg = PETSC_FALSE, set;
4682ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;
4782ea9b62SBarry Smith 
4882ea9b62SBarry Smith   PetscFunctionBegin;
4982ea9b62SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr);
5082ea9b62SBarry Smith   ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr);
519a0502c6SHåkon Strandenes   ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr);
52ee8b9732SVaclav Hapla   ierr = PetscOptionsBool("-viewer_hdf5_collective","Enable collective transfer mode","PetscViewerHDF5SetCollective",flg,&flg,&set);CHKERRQ(ierr);
53ee8b9732SVaclav Hapla   if (set) {ierr = PetscViewerHDF5SetCollective(v,flg);CHKERRQ(ierr);}
5482ea9b62SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
5582ea9b62SBarry Smith   PetscFunctionReturn(0);
5682ea9b62SBarry Smith }
5782ea9b62SBarry Smith 
581b793a25SVaclav Hapla static PetscErrorCode PetscViewerView_HDF5(PetscViewer v,PetscViewer viewer)
591b793a25SVaclav Hapla {
601b793a25SVaclav Hapla   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*)v->data;
6103fa8834SVaclav Hapla   PetscBool         flg;
621b793a25SVaclav Hapla   PetscErrorCode    ierr;
631b793a25SVaclav Hapla 
641b793a25SVaclav Hapla   PetscFunctionBegin;
651b793a25SVaclav Hapla   if (hdf5->filename) {
661b793a25SVaclav Hapla     ierr = PetscViewerASCIIPrintf(viewer,"Filename: %s\n",hdf5->filename);CHKERRQ(ierr);
671b793a25SVaclav Hapla   }
681b793a25SVaclav Hapla   ierr = PetscViewerASCIIPrintf(viewer,"Vectors with blocksize 1 saved as 2D datasets: %s\n",PetscBools[hdf5->basedimension2]);CHKERRQ(ierr);
691b793a25SVaclav Hapla   ierr = PetscViewerASCIIPrintf(viewer,"Enforce single precision storage: %s\n",PetscBools[hdf5->spoutput]);CHKERRQ(ierr);
7003fa8834SVaclav Hapla   ierr = PetscViewerHDF5GetCollective(v,&flg);CHKERRQ(ierr);
7103fa8834SVaclav Hapla   ierr = PetscViewerASCIIPrintf(viewer,"MPI-IO transfer mode: %s\n",flg ? "collective" : "independent");CHKERRQ(ierr);
721b793a25SVaclav Hapla   PetscFunctionReturn(0);
731b793a25SVaclav Hapla }
741b793a25SVaclav Hapla 
755c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
765c6c1daeSBarry Smith {
775c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
785c6c1daeSBarry Smith   PetscErrorCode   ierr;
795c6c1daeSBarry Smith 
805c6c1daeSBarry Smith   PetscFunctionBegin;
815c6c1daeSBarry Smith   ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);
82729ab6d9SBarry Smith   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
835c6c1daeSBarry Smith   PetscFunctionReturn(0);
845c6c1daeSBarry Smith }
855c6c1daeSBarry Smith 
869b1403beSVaclav Hapla static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
875c6c1daeSBarry Smith {
885c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
895c6c1daeSBarry Smith   PetscErrorCode   ierr;
905c6c1daeSBarry Smith 
915c6c1daeSBarry Smith   PetscFunctionBegin;
929c5072fbSVaclav Hapla   PetscStackCallHDF5(H5Pclose,(hdf5->dxpl_id));
935c6c1daeSBarry Smith   ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr);
945c6c1daeSBarry Smith   while (hdf5->groups) {
957d964842SVaclav Hapla     PetscViewerHDF5GroupList *tmp = hdf5->groups->next;
965c6c1daeSBarry Smith 
975c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups->name);CHKERRQ(ierr);
985c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups);CHKERRQ(ierr);
995c6c1daeSBarry Smith     hdf5->groups = tmp;
1005c6c1daeSBarry Smith   }
1015c6c1daeSBarry Smith   ierr = PetscFree(hdf5);CHKERRQ(ierr);
1020b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
103d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr);
1040b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
105058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr);
106058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr);
1075c6c1daeSBarry Smith   PetscFunctionReturn(0);
1085c6c1daeSBarry Smith }
1095c6c1daeSBarry Smith 
1109b1403beSVaclav Hapla static PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
1115c6c1daeSBarry Smith {
1125c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1135c6c1daeSBarry Smith 
1145c6c1daeSBarry Smith   PetscFunctionBegin;
1155c6c1daeSBarry Smith   hdf5->btype = type;
1165c6c1daeSBarry Smith   PetscFunctionReturn(0);
1175c6c1daeSBarry Smith }
1185c6c1daeSBarry Smith 
1190b62783dSVaclav Hapla static PetscErrorCode  PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
1200b62783dSVaclav Hapla {
1210b62783dSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1220b62783dSVaclav Hapla 
1230b62783dSVaclav Hapla   PetscFunctionBegin;
1240b62783dSVaclav Hapla   *type = hdf5->btype;
1250b62783dSVaclav Hapla   PetscFunctionReturn(0);
1260b62783dSVaclav Hapla }
1270b62783dSVaclav Hapla 
1289b1403beSVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
12982ea9b62SBarry Smith {
13082ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
13182ea9b62SBarry Smith 
13282ea9b62SBarry Smith   PetscFunctionBegin;
13382ea9b62SBarry Smith   hdf5->basedimension2 = flg;
13482ea9b62SBarry Smith   PetscFunctionReturn(0);
13582ea9b62SBarry Smith }
13682ea9b62SBarry Smith 
137df863907SAlex Fikl /*@
13882ea9b62SBarry Smith      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
13982ea9b62SBarry Smith        dimension of 2.
14082ea9b62SBarry Smith 
14182ea9b62SBarry Smith     Logically Collective on PetscViewer
14282ea9b62SBarry Smith 
14382ea9b62SBarry Smith   Input Parameters:
14482ea9b62SBarry Smith +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
14582ea9b62SBarry 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
14682ea9b62SBarry Smith 
14782ea9b62SBarry Smith   Options Database:
14882ea9b62SBarry 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
14982ea9b62SBarry Smith 
15095452b02SPatrick Sanan   Notes:
15195452b02SPatrick 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
15282ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
15382ea9b62SBarry Smith 
15482ea9b62SBarry Smith   Level: intermediate
15582ea9b62SBarry Smith 
15682ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
15782ea9b62SBarry Smith 
15882ea9b62SBarry Smith @*/
15982ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
16082ea9b62SBarry Smith {
16182ea9b62SBarry Smith   PetscErrorCode ierr;
16282ea9b62SBarry Smith 
16382ea9b62SBarry Smith   PetscFunctionBegin;
16482ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
16582ea9b62SBarry Smith   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
16682ea9b62SBarry Smith   PetscFunctionReturn(0);
16782ea9b62SBarry Smith }
16882ea9b62SBarry Smith 
169df863907SAlex Fikl /*@
17082ea9b62SBarry Smith      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
17182ea9b62SBarry Smith        dimension of 2.
17282ea9b62SBarry Smith 
17382ea9b62SBarry Smith     Logically Collective on PetscViewer
17482ea9b62SBarry Smith 
17582ea9b62SBarry Smith   Input Parameter:
17682ea9b62SBarry Smith .  viewer - the PetscViewer, must be of type HDF5
17782ea9b62SBarry Smith 
17882ea9b62SBarry Smith   Output Parameter:
17982ea9b62SBarry 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
18082ea9b62SBarry Smith 
18195452b02SPatrick Sanan   Notes:
18295452b02SPatrick 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
18382ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
18482ea9b62SBarry Smith 
18582ea9b62SBarry Smith   Level: intermediate
18682ea9b62SBarry Smith 
18782ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
18882ea9b62SBarry Smith 
18982ea9b62SBarry Smith @*/
19082ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
19182ea9b62SBarry Smith {
19282ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
19382ea9b62SBarry Smith 
19482ea9b62SBarry Smith   PetscFunctionBegin;
19582ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
19682ea9b62SBarry Smith   *flg = hdf5->basedimension2;
19782ea9b62SBarry Smith   PetscFunctionReturn(0);
19882ea9b62SBarry Smith }
19982ea9b62SBarry Smith 
2009b1403beSVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
2019a0502c6SHåkon Strandenes {
2029a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2039a0502c6SHåkon Strandenes 
2049a0502c6SHåkon Strandenes   PetscFunctionBegin;
2059a0502c6SHåkon Strandenes   hdf5->spoutput = flg;
2069a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2079a0502c6SHåkon Strandenes }
2089a0502c6SHåkon Strandenes 
209df863907SAlex Fikl /*@
2109a0502c6SHåkon Strandenes      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
2119a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2129a0502c6SHåkon Strandenes 
2139a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2149a0502c6SHåkon Strandenes 
2159a0502c6SHåkon Strandenes   Input Parameters:
2169a0502c6SHåkon Strandenes +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
2179a0502c6SHåkon Strandenes -  flg - if PETSC_TRUE the data will be written to disk with single precision
2189a0502c6SHåkon Strandenes 
2199a0502c6SHåkon Strandenes   Options Database:
2209a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
2219a0502c6SHåkon Strandenes 
22295452b02SPatrick Sanan   Notes:
22395452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2249a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2259a0502c6SHåkon Strandenes 
2269a0502c6SHåkon Strandenes   Level: intermediate
2279a0502c6SHåkon Strandenes 
2289a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2299a0502c6SHåkon Strandenes           PetscReal
2309a0502c6SHåkon Strandenes 
2319a0502c6SHåkon Strandenes @*/
2329a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
2339a0502c6SHåkon Strandenes {
2349a0502c6SHåkon Strandenes   PetscErrorCode ierr;
2359a0502c6SHåkon Strandenes 
2369a0502c6SHåkon Strandenes   PetscFunctionBegin;
2379a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2389a0502c6SHåkon Strandenes   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
2399a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2409a0502c6SHåkon Strandenes }
2419a0502c6SHåkon Strandenes 
242df863907SAlex Fikl /*@
2439a0502c6SHåkon Strandenes      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
2449a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2459a0502c6SHåkon Strandenes 
2469a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2479a0502c6SHåkon Strandenes 
2489a0502c6SHåkon Strandenes   Input Parameter:
2499a0502c6SHåkon Strandenes .  viewer - the PetscViewer, must be of type HDF5
2509a0502c6SHåkon Strandenes 
2519a0502c6SHåkon Strandenes   Output Parameter:
2529a0502c6SHåkon Strandenes .  flg - if PETSC_TRUE the data will be written to disk with single precision
2539a0502c6SHåkon Strandenes 
25495452b02SPatrick Sanan   Notes:
25595452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2569a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2579a0502c6SHåkon Strandenes 
2589a0502c6SHåkon Strandenes   Level: intermediate
2599a0502c6SHåkon Strandenes 
2609a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2619a0502c6SHåkon Strandenes           PetscReal
2629a0502c6SHåkon Strandenes 
2639a0502c6SHåkon Strandenes @*/
2649a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
2659a0502c6SHåkon Strandenes {
2669a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2679a0502c6SHåkon Strandenes 
2689a0502c6SHåkon Strandenes   PetscFunctionBegin;
2699a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2709a0502c6SHåkon Strandenes   *flg = hdf5->spoutput;
2719a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2729a0502c6SHåkon Strandenes }
2739a0502c6SHåkon Strandenes 
274ee8b9732SVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg)
275ee8b9732SVaclav Hapla {
276ee8b9732SVaclav Hapla   PetscFunctionBegin;
277ee8b9732SVaclav Hapla   /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions
278ee8b9732SVaclav Hapla      - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */
279c796909eSBarry Smith #if H5_VERSION_GE(1,10,3) && defined(H5_HAVE_PARALLEL)
280ee8b9732SVaclav Hapla   {
281ee8b9732SVaclav Hapla     PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
282ee8b9732SVaclav Hapla     PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT));
283ee8b9732SVaclav Hapla   }
284ee8b9732SVaclav Hapla #else
285ee8b9732SVaclav Hapla   if (flg) {
286ee8b9732SVaclav Hapla     PetscErrorCode ierr;
287c796909eSBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)viewer), "Warning: PetscViewerHDF5SetCollective(viewer,PETSC_TRUE) is ignored for HDF5 versions prior to 1.10.3 or if built without MPI support\n");CHKERRQ(ierr);
288ee8b9732SVaclav Hapla   }
289ee8b9732SVaclav Hapla #endif
290ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
291ee8b9732SVaclav Hapla }
292ee8b9732SVaclav Hapla 
293ee8b9732SVaclav Hapla /*@
294ee8b9732SVaclav Hapla   PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes.
295ee8b9732SVaclav Hapla 
296ee8b9732SVaclav Hapla   Logically Collective; flg must contain common value
297ee8b9732SVaclav Hapla 
298ee8b9732SVaclav Hapla   Input Parameters:
299ee8b9732SVaclav Hapla + viewer - the PetscViewer; if it is not hdf5 then this command is ignored
300ee8b9732SVaclav Hapla - flg - PETSC_TRUE for collective mode; PETSC_FALSE for independent mode (default)
301ee8b9732SVaclav Hapla 
302ee8b9732SVaclav Hapla   Options Database:
303ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers
304ee8b9732SVaclav Hapla 
305ee8b9732SVaclav Hapla   Notes:
306ee8b9732SVaclav Hapla   Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better.
30753effed7SBarry Smith   However, this works correctly only since HDF5 1.10.3 and if HDF5 is installed for MPI; hence, we ignore this setting for older versions.
308ee8b9732SVaclav Hapla 
309ee8b9732SVaclav Hapla   Developer notes:
310ee8b9732SVaclav Hapla   In the HDF5 layer, PETSC_TRUE / PETSC_FALSE means H5Pset_dxpl_mpio() is called with H5FD_MPIO_COLLECTIVE / H5FD_MPIO_INDEPENDENT, respectively.
311ee8b9732SVaclav Hapla   This in turn means use of MPI_File_{read,write}_all /  MPI_File_{read,write} in the MPI-IO layer, respectively.
312ee8b9732SVaclav Hapla   See HDF5 documentation and MPI-IO documentation for details.
313ee8b9732SVaclav Hapla 
314ee8b9732SVaclav Hapla   Level: intermediate
315ee8b9732SVaclav Hapla 
316ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5GetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open()
317ee8b9732SVaclav Hapla 
318ee8b9732SVaclav Hapla @*/
319ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg)
320ee8b9732SVaclav Hapla {
321ee8b9732SVaclav Hapla   PetscErrorCode ierr;
322ee8b9732SVaclav Hapla 
323ee8b9732SVaclav Hapla   PetscFunctionBegin;
324ee8b9732SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
325ee8b9732SVaclav Hapla   PetscValidLogicalCollectiveBool(viewer,flg,2);
326ee8b9732SVaclav Hapla   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
327ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
328ee8b9732SVaclav Hapla }
329ee8b9732SVaclav Hapla 
330ee8b9732SVaclav Hapla static PetscErrorCode  PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg)
331ee8b9732SVaclav Hapla {
332c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL)
333ee8b9732SVaclav Hapla   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*) viewer->data;
334ee8b9732SVaclav Hapla   H5FD_mpio_xfer_t  mode;
335c796909eSBarry Smith #endif
336ee8b9732SVaclav Hapla 
337ee8b9732SVaclav Hapla   PetscFunctionBegin;
338c796909eSBarry Smith #if !defined(H5_HAVE_PARALLEL)
339c796909eSBarry Smith   *flg = PETSC_FALSE;
340c796909eSBarry Smith #else
341ee8b9732SVaclav Hapla   PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode));
342ee8b9732SVaclav Hapla   *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE;
343c796909eSBarry Smith #endif
344ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
345ee8b9732SVaclav Hapla }
346ee8b9732SVaclav Hapla 
347ee8b9732SVaclav Hapla /*@
348ee8b9732SVaclav Hapla   PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes.
349ee8b9732SVaclav Hapla 
350ee8b9732SVaclav Hapla   Not Collective
351ee8b9732SVaclav Hapla 
352ee8b9732SVaclav Hapla   Input Parameters:
353ee8b9732SVaclav Hapla . viewer - the HDF5 PetscViewer
354ee8b9732SVaclav Hapla 
355ee8b9732SVaclav Hapla   Output Parameters:
356ee8b9732SVaclav Hapla . flg - the flag
357ee8b9732SVaclav Hapla 
358ee8b9732SVaclav Hapla   Level: intermediate
359ee8b9732SVaclav Hapla 
360ee8b9732SVaclav Hapla   Notes:
361c796909eSBarry Smith   This setting works correctly only since HDF5 1.10.3 and if HDF5 was installed for MPI. For older versions, PETSC_FALSE will be always returned.
362ee8b9732SVaclav Hapla   For more details, see PetscViewerHDF5SetCollective().
363ee8b9732SVaclav Hapla 
364ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open()
365ee8b9732SVaclav Hapla 
366ee8b9732SVaclav Hapla @*/
367ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg)
368ee8b9732SVaclav Hapla {
369ee8b9732SVaclav Hapla   PetscErrorCode ierr;
370ee8b9732SVaclav Hapla 
371ee8b9732SVaclav Hapla   PetscFunctionBegin;
372ee8b9732SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
373534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg,2);
374ee8b9732SVaclav Hapla 
375ee8b9732SVaclav Hapla   ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));CHKERRQ(ierr);
376ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
377ee8b9732SVaclav Hapla }
378ee8b9732SVaclav Hapla 
3799b1403beSVaclav Hapla static PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
3805c6c1daeSBarry Smith {
3815c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3825c6c1daeSBarry Smith   hid_t             plist_id;
3835c6c1daeSBarry Smith   PetscErrorCode    ierr;
3845c6c1daeSBarry Smith 
3855c6c1daeSBarry Smith   PetscFunctionBegin;
386f683cc58SToby Isaac   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
387f683cc58SToby Isaac   if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);}
3885c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
3895c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
390729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
391c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL)
392d28bfa55SVaclav Hapla   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL));
393c796909eSBarry Smith #endif
3945c6c1daeSBarry Smith   /* Create or open the file collectively */
3955c6c1daeSBarry Smith   switch (hdf5->btype) {
3965c6c1daeSBarry Smith   case FILE_MODE_READ:
397729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
3985c6c1daeSBarry Smith     break;
3995c6c1daeSBarry Smith   case FILE_MODE_APPEND:
4007e4fd573SVaclav Hapla   case FILE_MODE_UPDATE:
4017e4fd573SVaclav Hapla   {
4027e4fd573SVaclav Hapla     PetscBool flg;
4037e4fd573SVaclav Hapla     ierr = PetscTestFile(hdf5->filename, 'r', &flg);CHKERRQ(ierr);
4047e4fd573SVaclav Hapla     if (flg) PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
4057e4fd573SVaclav Hapla     else     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id));
4065c6c1daeSBarry Smith     break;
4077e4fd573SVaclav Hapla   }
4085c6c1daeSBarry Smith   case FILE_MODE_WRITE:
409729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
4105c6c1daeSBarry Smith     break;
4117e4fd573SVaclav Hapla   case FILE_MODE_UNDEFINED:
4127e4fd573SVaclav Hapla     SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
4135c6c1daeSBarry Smith   default:
4147e4fd573SVaclav Hapla     SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP, "Unsupported file mode %s",PetscFileModes[hdf5->btype]);
4155c6c1daeSBarry Smith   }
4165c6c1daeSBarry Smith   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
417729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
4185c6c1daeSBarry Smith   PetscFunctionReturn(0);
4195c6c1daeSBarry Smith }
4205c6c1daeSBarry Smith 
421d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
422d1232d7fSToby Isaac {
423d1232d7fSToby Isaac   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
424d1232d7fSToby Isaac 
425d1232d7fSToby Isaac   PetscFunctionBegin;
426d1232d7fSToby Isaac   *name = vhdf5->filename;
427d1232d7fSToby Isaac   PetscFunctionReturn(0);
428d1232d7fSToby Isaac }
429d1232d7fSToby Isaac 
430b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
431b723ab35SVaclav Hapla {
4325dc64a97SVaclav Hapla   /*
433b723ab35SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
434b723ab35SVaclav Hapla   PetscErrorCode   ierr;
4355dc64a97SVaclav Hapla   */
436b723ab35SVaclav Hapla 
437b723ab35SVaclav Hapla   PetscFunctionBegin;
438b723ab35SVaclav Hapla   PetscFunctionReturn(0);
439b723ab35SVaclav Hapla }
440b723ab35SVaclav Hapla 
4418556b5ebSBarry Smith /*MC
4428556b5ebSBarry Smith    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
4438556b5ebSBarry Smith 
4448556b5ebSBarry Smith .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
4458556b5ebSBarry Smith            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
4468556b5ebSBarry Smith            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
4478556b5ebSBarry Smith            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
4488556b5ebSBarry Smith 
4491b266c99SBarry Smith   Level: beginner
4508556b5ebSBarry Smith M*/
451d1232d7fSToby Isaac 
4528cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
4535c6c1daeSBarry Smith {
4545c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
4555c6c1daeSBarry Smith   PetscErrorCode   ierr;
4565c6c1daeSBarry Smith 
4575c6c1daeSBarry Smith   PetscFunctionBegin;
458b00a9115SJed Brown   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
4595c6c1daeSBarry Smith 
4605c6c1daeSBarry Smith   v->data                = (void*) hdf5;
4615c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
46282ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
463b723ab35SVaclav Hapla   v->ops->setup          = PetscViewerSetUp_HDF5;
4641b793a25SVaclav Hapla   v->ops->view           = PetscViewerView_HDF5;
465908793a3SLisandro Dalcin   v->ops->flush          = NULL;
4667e4fd573SVaclav Hapla   hdf5->btype            = FILE_MODE_UNDEFINED;
467908793a3SLisandro Dalcin   hdf5->filename         = NULL;
4685c6c1daeSBarry Smith   hdf5->timestep         = -1;
4690298fd71SBarry Smith   hdf5->groups           = NULL;
4705c6c1daeSBarry Smith 
4719c5072fbSVaclav Hapla   PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER));
4729c5072fbSVaclav Hapla 
4730b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
474d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr);
4750b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
4760b62783dSVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr);
47782ea9b62SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
4789a0502c6SHåkon Strandenes   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr);
479ee8b9732SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);CHKERRQ(ierr);
480ee8b9732SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);CHKERRQ(ierr);
4815c6c1daeSBarry Smith   PetscFunctionReturn(0);
4825c6c1daeSBarry Smith }
4835c6c1daeSBarry Smith 
4845c6c1daeSBarry Smith /*@C
4855c6c1daeSBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
4865c6c1daeSBarry Smith 
487d083f849SBarry Smith    Collective
4885c6c1daeSBarry Smith 
4895c6c1daeSBarry Smith    Input Parameters:
4905c6c1daeSBarry Smith +  comm - MPI communicator
4915c6c1daeSBarry Smith .  name - name of file
4925c6c1daeSBarry Smith -  type - type of file
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:
498a2b725a8SWilliam Gropp +  -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
499a2b725a8SWilliam Gropp -  -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 
5037e4fd573SVaclav Hapla    Notes:
5047e4fd573SVaclav Hapla    Reading is always available, regardless of the mode. Available modes are
5057e4fd573SVaclav Hapla +  FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY]
5067e4fd573SVaclav Hapla .  FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC]
5077e4fd573SVaclav Hapla .  FILE_MODE_APPEND - if file exists, keep existing contents [H5Fopen() with H5F_ACC_RDWR], else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_EXCL]
5087e4fd573SVaclav Hapla -  FILE_MODE_UPDATE - same as FILE_MODE_APPEND
5097e4fd573SVaclav Hapla 
5107e4fd573SVaclav Hapla    In case of FILE_MODE_APPEND / FILE_MODE_UPDATE, any stored object (dataset, attribute) can be selectively ovewritten if the same fully qualified name (/group/path/to/object) is specified.
5117e4fd573SVaclav Hapla 
5125c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
5135c6c1daeSBarry Smith 
5146a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
5159a0502c6SHåkon Strandenes           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
516a56f64adSBarry Smith           MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
5175c6c1daeSBarry Smith @*/
5185c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
5195c6c1daeSBarry Smith {
5205c6c1daeSBarry Smith   PetscErrorCode ierr;
5215c6c1daeSBarry Smith 
5225c6c1daeSBarry Smith   PetscFunctionBegin;
5235c6c1daeSBarry Smith   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
5245c6c1daeSBarry Smith   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
5255c6c1daeSBarry Smith   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
5265c6c1daeSBarry Smith   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
527b663d2d2SVaclav Hapla   ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr);
5285c6c1daeSBarry Smith   PetscFunctionReturn(0);
5295c6c1daeSBarry Smith }
5305c6c1daeSBarry Smith 
5315c6c1daeSBarry Smith /*@C
5325c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
5335c6c1daeSBarry Smith 
5345c6c1daeSBarry Smith   Not collective
5355c6c1daeSBarry Smith 
5365c6c1daeSBarry Smith   Input Parameter:
5375c6c1daeSBarry Smith . viewer - the PetscViewer
5385c6c1daeSBarry Smith 
5395c6c1daeSBarry Smith   Output Parameter:
5405c6c1daeSBarry Smith . file_id - The file id
5415c6c1daeSBarry Smith 
5425c6c1daeSBarry Smith   Level: intermediate
5435c6c1daeSBarry Smith 
5445c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open()
5455c6c1daeSBarry Smith @*/
5465c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
5475c6c1daeSBarry Smith {
5485c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5495c6c1daeSBarry Smith 
5505c6c1daeSBarry Smith   PetscFunctionBegin;
5515c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5525c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
5535c6c1daeSBarry Smith   PetscFunctionReturn(0);
5545c6c1daeSBarry Smith }
5555c6c1daeSBarry Smith 
5565c6c1daeSBarry Smith /*@C
5575c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
5585c6c1daeSBarry Smith 
5595c6c1daeSBarry Smith   Not collective
5605c6c1daeSBarry Smith 
5615c6c1daeSBarry Smith   Input Parameters:
5625c6c1daeSBarry Smith + viewer - the PetscViewer
5635c6c1daeSBarry Smith - name - The group name
5645c6c1daeSBarry Smith 
5655c6c1daeSBarry Smith   Level: intermediate
5665c6c1daeSBarry Smith 
5672e361344SVaclav Hapla   Notes:
5682e361344SVaclav Hapla   This is designed to mnemonically resemble the Unix cd command.
5692e361344SVaclav Hapla   + If name begins with '/', it is interpreted as an absolute path fully replacing current group, otherwise it is taken as relative to the current group.
5702e361344SVaclav Hapla   . NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/".
5712e361344SVaclav Hapla   - "." means the current group is pushed again.
5722e361344SVaclav Hapla 
5732e361344SVaclav Hapla   Example:
5742e361344SVaclav Hapla   Suppose the current group is "/a".
5752e361344SVaclav Hapla   + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/".
5762e361344SVaclav Hapla   . If name is ".", then the new group will be "/a".
5772e361344SVaclav Hapla   . If name is "b", then the new group will be "/a/b".
5782e361344SVaclav Hapla   - If name is "/b", then the new group will be "/b".
5792e361344SVaclav Hapla 
5802e361344SVaclav Hapla   Developer Notes:
5812e361344SVaclav Hapla   The root group "/" is internally stored as NULL.
582820fc2d1SVaclav Hapla 
583874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
5845c6c1daeSBarry Smith @*/
585be7aa430SVaclav Hapla PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[])
5865c6c1daeSBarry Smith {
5875c6c1daeSBarry Smith   PetscViewer_HDF5          *hdf5 = (PetscViewer_HDF5*) viewer->data;
5887d964842SVaclav Hapla   PetscViewerHDF5GroupList  *groupNode;
5892e361344SVaclav Hapla   size_t                    i,len;
5902e361344SVaclav Hapla   char                      buf[PETSC_MAX_PATH_LEN];
5912e361344SVaclav Hapla   const char                *gname;
5925c6c1daeSBarry Smith   PetscErrorCode            ierr;
5935c6c1daeSBarry Smith 
5945c6c1daeSBarry Smith   PetscFunctionBegin;
5955c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
596820fc2d1SVaclav Hapla   if (name) PetscValidCharPointer(name,2);
597820fc2d1SVaclav Hapla   ierr = PetscStrlen(name, &len);CHKERRQ(ierr);
5982e361344SVaclav Hapla   gname = NULL;
5992e361344SVaclav Hapla   if (len) {
6002e361344SVaclav Hapla     if (len == 1 && name[0] == '.') {
6012e361344SVaclav Hapla       /* use current name */
6022e361344SVaclav Hapla       gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL;
6032e361344SVaclav Hapla     } else if (name[0] == '/') {
6042e361344SVaclav Hapla       /* absolute */
6052e361344SVaclav Hapla       for (i=1; i<len; i++) {
6062e361344SVaclav Hapla         if (name[i] != '/') {
6072e361344SVaclav Hapla           gname = name;
6082e361344SVaclav Hapla           break;
6092e361344SVaclav Hapla         }
6102e361344SVaclav Hapla       }
6112e361344SVaclav Hapla     } else {
6122e361344SVaclav Hapla       /* relative */
6132e361344SVaclav Hapla       const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : "";
6142e361344SVaclav Hapla       ierr = PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name);CHKERRQ(ierr);
6152e361344SVaclav Hapla       gname = buf;
6162e361344SVaclav Hapla     }
6172e361344SVaclav Hapla   }
61895dccacaSBarry Smith   ierr = PetscNew(&groupNode);CHKERRQ(ierr);
6192e361344SVaclav Hapla   ierr = PetscStrallocpy(gname, (char**) &groupNode->name);CHKERRQ(ierr);
6205c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
6215c6c1daeSBarry Smith   hdf5->groups    = groupNode;
6225c6c1daeSBarry Smith   PetscFunctionReturn(0);
6235c6c1daeSBarry Smith }
6245c6c1daeSBarry Smith 
6253ef9c667SSatish Balay /*@
6265c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
6275c6c1daeSBarry Smith 
6285c6c1daeSBarry Smith   Not collective
6295c6c1daeSBarry Smith 
6305c6c1daeSBarry Smith   Input Parameter:
6315c6c1daeSBarry Smith . viewer - the PetscViewer
6325c6c1daeSBarry Smith 
6335c6c1daeSBarry Smith   Level: intermediate
6345c6c1daeSBarry Smith 
635874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
6365c6c1daeSBarry Smith @*/
6375c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
6385c6c1daeSBarry Smith {
6395c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
6407d964842SVaclav Hapla   PetscViewerHDF5GroupList *groupNode;
6415c6c1daeSBarry Smith   PetscErrorCode   ierr;
6425c6c1daeSBarry Smith 
6435c6c1daeSBarry Smith   PetscFunctionBegin;
6445c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
64582f516ccSBarry Smith   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
6465c6c1daeSBarry Smith   groupNode    = hdf5->groups;
6475c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
6485c6c1daeSBarry Smith   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
6495c6c1daeSBarry Smith   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
6505c6c1daeSBarry Smith   PetscFunctionReturn(0);
6515c6c1daeSBarry Smith }
6525c6c1daeSBarry Smith 
6535c6c1daeSBarry Smith /*@C
654874270d9SVaclav Hapla   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup().
655874270d9SVaclav Hapla   If none has been assigned, returns NULL.
6565c6c1daeSBarry Smith 
6575c6c1daeSBarry Smith   Not collective
6585c6c1daeSBarry Smith 
6595c6c1daeSBarry Smith   Input Parameter:
6605c6c1daeSBarry Smith . viewer - the PetscViewer
6615c6c1daeSBarry Smith 
6625c6c1daeSBarry Smith   Output Parameter:
6635c6c1daeSBarry Smith . name - The group name
6645c6c1daeSBarry Smith 
6655c6c1daeSBarry Smith   Level: intermediate
6665c6c1daeSBarry Smith 
667874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
6685c6c1daeSBarry Smith @*/
669be7aa430SVaclav Hapla PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[])
6705c6c1daeSBarry Smith {
6715c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
6725c6c1daeSBarry Smith 
6735c6c1daeSBarry Smith   PetscFunctionBegin;
6745c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
675c959eef4SJed Brown   PetscValidPointer(name,2);
676a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
6770298fd71SBarry Smith   else *name = NULL;
6785c6c1daeSBarry Smith   PetscFunctionReturn(0);
6795c6c1daeSBarry Smith }
6805c6c1daeSBarry Smith 
6818c557b5aSVaclav Hapla /*@
682874270d9SVaclav Hapla   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(),
683874270d9SVaclav Hapla   and return this group's ID and file ID.
684874270d9SVaclav Hapla   If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID.
685874270d9SVaclav Hapla 
686874270d9SVaclav Hapla   Not collective
687874270d9SVaclav Hapla 
688874270d9SVaclav Hapla   Input Parameter:
689874270d9SVaclav Hapla . viewer - the PetscViewer
690874270d9SVaclav Hapla 
691874270d9SVaclav Hapla   Output Parameter:
692874270d9SVaclav Hapla + fileId - The HDF5 file ID
693874270d9SVaclav Hapla - groupId - The HDF5 group ID
694874270d9SVaclav Hapla 
695e74616beSVaclav Hapla   Notes:
696e74616beSVaclav Hapla   If the viewer is writable, the group is created if it doesn't exist yet.
697e74616beSVaclav Hapla 
698874270d9SVaclav Hapla   Level: intermediate
699874270d9SVaclav Hapla 
700874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
701874270d9SVaclav Hapla @*/
70254dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
70354dbf706SVaclav Hapla {
70490e03537SVaclav Hapla   hid_t          file_id;
70590e03537SVaclav Hapla   H5O_type_t     type;
70654dbf706SVaclav Hapla   const char     *groupName = NULL;
707e74616beSVaclav Hapla   PetscBool      create;
70854dbf706SVaclav Hapla   PetscErrorCode ierr;
70954dbf706SVaclav Hapla 
71054dbf706SVaclav Hapla   PetscFunctionBegin;
711e74616beSVaclav Hapla   ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr);
71254dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr);
71354dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr);
714e74616beSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr);
71590e03537SVaclav 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);
71690e03537SVaclav Hapla   PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT));
71754dbf706SVaclav Hapla   *fileId  = file_id;
71854dbf706SVaclav Hapla   PetscFunctionReturn(0);
71954dbf706SVaclav Hapla }
72054dbf706SVaclav Hapla 
7213ef9c667SSatish Balay /*@
722d7dd068bSVaclav Hapla   PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing.
7235c6c1daeSBarry Smith 
7245c6c1daeSBarry Smith   Not collective
7255c6c1daeSBarry Smith 
7265c6c1daeSBarry Smith   Input Parameter:
7275c6c1daeSBarry Smith . viewer - the PetscViewer
7285c6c1daeSBarry Smith 
7295c6c1daeSBarry Smith   Level: intermediate
7305c6c1daeSBarry Smith 
731d7dd068bSVaclav Hapla   Notes:
732d7dd068bSVaclav Hapla   On first PetscViewerHDF5PushTimestepping(), the initial time step is set to 0.
733d7dd068bSVaclav Hapla   Next timesteps can then be set using PetscViewerHDF5IncrementTimestep() or PetscViewerHDF5SetTimestep().
734d7dd068bSVaclav Hapla   Current timestep value determines which timestep is read from or written to any dataset on the next HDF5 I/O operation [e.g. VecView()].
735d7dd068bSVaclav Hapla   Use PetscViewerHDF5PopTimestepping() to deactivate timestepping mode; calling it by the end of the program is NOT mandatory.
736d7dd068bSVaclav Hapla   Current timestep is remembered between PetscViewerHDF5PopTimestepping() and the next PetscViewerHDF5PushTimestepping().
737d7dd068bSVaclav Hapla 
738d7dd068bSVaclav Hapla   If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again.
739d7dd068bSVaclav Hapla   Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error.
740d7dd068bSVaclav Hapla 
741d7dd068bSVaclav Hapla   Developer notes:
742d7dd068bSVaclav Hapla   Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true.
743d7dd068bSVaclav Hapla 
744d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PopTimestepping(), PetscViewerHDF5IsTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
745d7dd068bSVaclav Hapla @*/
746d7dd068bSVaclav Hapla PetscErrorCode  PetscViewerHDF5PushTimestepping(PetscViewer viewer)
747d7dd068bSVaclav Hapla {
748d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
749d7dd068bSVaclav Hapla 
750d7dd068bSVaclav Hapla   PetscFunctionBegin;
751d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
752d7dd068bSVaclav Hapla   if (hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed");
753d7dd068bSVaclav Hapla   hdf5->timestepping = PETSC_TRUE;
754d7dd068bSVaclav Hapla   if (hdf5->timestep < 0) hdf5->timestep = 0;
755d7dd068bSVaclav Hapla   PetscFunctionReturn(0);
756d7dd068bSVaclav Hapla }
757d7dd068bSVaclav Hapla 
758d7dd068bSVaclav Hapla /*@
759d7dd068bSVaclav Hapla   PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing.
760d7dd068bSVaclav Hapla 
761d7dd068bSVaclav Hapla   Not collective
762d7dd068bSVaclav Hapla 
763d7dd068bSVaclav Hapla   Input Parameter:
764d7dd068bSVaclav Hapla . viewer - the PetscViewer
765d7dd068bSVaclav Hapla 
766d7dd068bSVaclav Hapla   Level: intermediate
767d7dd068bSVaclav Hapla 
768d7dd068bSVaclav Hapla   Notes:
769d7dd068bSVaclav Hapla   See PetscViewerHDF5PushTimestepping() for details.
770d7dd068bSVaclav Hapla 
771d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IsTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
772d7dd068bSVaclav Hapla @*/
773d7dd068bSVaclav Hapla PetscErrorCode  PetscViewerHDF5PopTimestepping(PetscViewer viewer)
774d7dd068bSVaclav Hapla {
775d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
776d7dd068bSVaclav Hapla 
777d7dd068bSVaclav Hapla   PetscFunctionBegin;
778d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
779d7dd068bSVaclav Hapla   if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
780d7dd068bSVaclav Hapla   hdf5->timestepping = PETSC_FALSE;
781d7dd068bSVaclav Hapla   PetscFunctionReturn(0);
782d7dd068bSVaclav Hapla }
783d7dd068bSVaclav Hapla 
784d7dd068bSVaclav Hapla /*@
785d7dd068bSVaclav Hapla   PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently.
786d7dd068bSVaclav Hapla 
787d7dd068bSVaclav Hapla   Not collective
788d7dd068bSVaclav Hapla 
789d7dd068bSVaclav Hapla   Input Parameter:
790d7dd068bSVaclav Hapla . viewer - the PetscViewer
791d7dd068bSVaclav Hapla 
792d7dd068bSVaclav Hapla   Output Parameter:
793d7dd068bSVaclav Hapla . flg - is timestepping active?
794d7dd068bSVaclav Hapla 
795d7dd068bSVaclav Hapla   Level: intermediate
796d7dd068bSVaclav Hapla 
797d7dd068bSVaclav Hapla   Notes:
798d7dd068bSVaclav Hapla   See PetscViewerHDF5PushTimestepping() for details.
799d7dd068bSVaclav Hapla 
800d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5PopTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
801d7dd068bSVaclav Hapla @*/
802d7dd068bSVaclav Hapla PetscErrorCode  PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg)
803d7dd068bSVaclav Hapla {
804d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
805d7dd068bSVaclav Hapla 
806d7dd068bSVaclav Hapla   PetscFunctionBegin;
807d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
808d7dd068bSVaclav Hapla   *flg = hdf5->timestepping;
809d7dd068bSVaclav Hapla   PetscFunctionReturn(0);
810d7dd068bSVaclav Hapla }
811d7dd068bSVaclav Hapla 
812d7dd068bSVaclav Hapla /*@
813d7dd068bSVaclav Hapla   PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time.
814d7dd068bSVaclav Hapla 
815d7dd068bSVaclav Hapla   Not collective
816d7dd068bSVaclav Hapla 
817d7dd068bSVaclav Hapla   Input Parameter:
818d7dd068bSVaclav Hapla . viewer - the PetscViewer
819d7dd068bSVaclav Hapla 
820d7dd068bSVaclav Hapla   Level: intermediate
821d7dd068bSVaclav Hapla 
822d7dd068bSVaclav Hapla   Notes:
823d7dd068bSVaclav Hapla   This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details.
824d7dd068bSVaclav Hapla 
825d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
8265c6c1daeSBarry Smith @*/
8275c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
8285c6c1daeSBarry Smith {
8295c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
8305c6c1daeSBarry Smith 
8315c6c1daeSBarry Smith   PetscFunctionBegin;
8325c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
833d7dd068bSVaclav Hapla   if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
8345c6c1daeSBarry Smith   ++hdf5->timestep;
8355c6c1daeSBarry Smith   PetscFunctionReturn(0);
8365c6c1daeSBarry Smith }
8375c6c1daeSBarry Smith 
8383ef9c667SSatish Balay /*@
839d7dd068bSVaclav Hapla   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time.
8405c6c1daeSBarry Smith 
841d7dd068bSVaclav Hapla   Logically collective
8425c6c1daeSBarry Smith 
8435c6c1daeSBarry Smith   Input Parameters:
8445c6c1daeSBarry Smith + viewer - the PetscViewer
845d7dd068bSVaclav Hapla - timestep - The timestep
8465c6c1daeSBarry Smith 
8475c6c1daeSBarry Smith   Level: intermediate
8485c6c1daeSBarry Smith 
849d7dd068bSVaclav Hapla   Notes:
850d7dd068bSVaclav Hapla   This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details.
851d7dd068bSVaclav Hapla 
852d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
8535c6c1daeSBarry Smith @*/
8545c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
8555c6c1daeSBarry Smith {
8565c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
8575c6c1daeSBarry Smith 
8585c6c1daeSBarry Smith   PetscFunctionBegin;
8595c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
860d7dd068bSVaclav Hapla   PetscValidLogicalCollectiveInt(viewer, timestep, 2);
861d7dd068bSVaclav Hapla   if (timestep < 0) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %D is negative", timestep);
862d7dd068bSVaclav Hapla   if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
8635c6c1daeSBarry Smith   hdf5->timestep = timestep;
8645c6c1daeSBarry Smith   PetscFunctionReturn(0);
8655c6c1daeSBarry Smith }
8665c6c1daeSBarry Smith 
8673ef9c667SSatish Balay /*@
8685c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
8695c6c1daeSBarry Smith 
8705c6c1daeSBarry Smith   Not collective
8715c6c1daeSBarry Smith 
8725c6c1daeSBarry Smith   Input Parameter:
8735c6c1daeSBarry Smith . viewer - the PetscViewer
8745c6c1daeSBarry Smith 
8755c6c1daeSBarry Smith   Output Parameter:
876d7dd068bSVaclav Hapla . timestep - The timestep
8775c6c1daeSBarry Smith 
8785c6c1daeSBarry Smith   Level: intermediate
8795c6c1daeSBarry Smith 
880d7dd068bSVaclav Hapla   Notes:
881d7dd068bSVaclav Hapla   This can be called only if the viewer is in the timestepping mode. See PetscViewerHDF5PushTimestepping() for details.
882d7dd068bSVaclav Hapla 
883d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
8845c6c1daeSBarry Smith @*/
8855c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
8865c6c1daeSBarry Smith {
8875c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
8885c6c1daeSBarry Smith 
8895c6c1daeSBarry Smith   PetscFunctionBegin;
8905c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
891d7dd068bSVaclav Hapla   PetscValidIntPointer(timestep,2);
892d7dd068bSVaclav Hapla   if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
8935c6c1daeSBarry Smith   *timestep = hdf5->timestep;
8945c6c1daeSBarry Smith   PetscFunctionReturn(0);
8955c6c1daeSBarry Smith }
8965c6c1daeSBarry Smith 
89736bce6e8SMatthew G. Knepley /*@C
89836bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
89936bce6e8SMatthew G. Knepley 
90036bce6e8SMatthew G. Knepley   Not collective
90136bce6e8SMatthew G. Knepley 
90236bce6e8SMatthew G. Knepley   Input Parameter:
90336bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
90436bce6e8SMatthew G. Knepley 
90536bce6e8SMatthew G. Knepley   Output Parameter:
90636bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
90736bce6e8SMatthew G. Knepley 
90836bce6e8SMatthew G. Knepley   Level: advanced
90936bce6e8SMatthew G. Knepley 
91036bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
91136bce6e8SMatthew G. Knepley @*/
91236bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
91336bce6e8SMatthew G. Knepley {
91436bce6e8SMatthew G. Knepley   PetscFunctionBegin;
91536bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
91636bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
91736bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_LLONG;
91836bce6e8SMatthew G. Knepley #else
91936bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_INT;
92036bce6e8SMatthew G. Knepley #endif
92136bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
92236bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
92336bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
924c9450970SBarry Smith   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_INT;
925de3ba262SVaclav Hapla   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_INT;
92636bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
92736bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
92836bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
9297e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
93036bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
93136bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
93236bce6e8SMatthew G. Knepley }
93336bce6e8SMatthew G. Knepley 
93436bce6e8SMatthew G. Knepley /*@C
93536bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
93636bce6e8SMatthew G. Knepley 
93736bce6e8SMatthew G. Knepley   Not collective
93836bce6e8SMatthew G. Knepley 
93936bce6e8SMatthew G. Knepley   Input Parameter:
94036bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
94136bce6e8SMatthew G. Knepley 
94236bce6e8SMatthew G. Knepley   Output Parameter:
94336bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
94436bce6e8SMatthew G. Knepley 
94536bce6e8SMatthew G. Knepley   Level: advanced
94636bce6e8SMatthew G. Knepley 
94736bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
94836bce6e8SMatthew G. Knepley @*/
94936bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
95036bce6e8SMatthew G. Knepley {
95136bce6e8SMatthew G. Knepley   PetscFunctionBegin;
95236bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
95336bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
95436bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
95536bce6e8SMatthew G. Knepley #else
95636bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
95736bce6e8SMatthew G. Knepley #endif
95836bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
95936bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
96036bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
96136bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
96236bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
96336bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
9647e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
96536bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
96636bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
96736bce6e8SMatthew G. Knepley }
96836bce6e8SMatthew G. Knepley 
969df863907SAlex Fikl /*@C
970b8ee85a1SVaclav Hapla  PetscViewerHDF5WriteAttribute - Write an attribute
97136bce6e8SMatthew G. Knepley 
972343a20b2SVaclav Hapla   Collective
973343a20b2SVaclav Hapla 
97436bce6e8SMatthew G. Knepley   Input Parameters:
97536bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer
9764302210dSVaclav Hapla . parent - The parent dataset/group name
97736bce6e8SMatthew G. Knepley . name   - The attribute name
97836bce6e8SMatthew G. Knepley . datatype - The attribute type
97936bce6e8SMatthew G. Knepley - value    - The attribute value
98036bce6e8SMatthew G. Knepley 
98136bce6e8SMatthew G. Knepley   Level: advanced
98236bce6e8SMatthew G. Knepley 
9834302210dSVaclav Hapla   Notes:
984343a20b2SVaclav Hapla   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group. NULL means the current pushed group.
9854302210dSVaclav Hapla 
986577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
98736bce6e8SMatthew G. Knepley @*/
9884302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
98936bce6e8SMatthew G. Knepley {
9904302210dSVaclav Hapla   char           *parentAbsPath;
99160568592SMatthew G. Knepley   hid_t          h5, dataspace, obj, attribute, dtype;
992080f144cSVaclav Hapla   PetscBool      has;
99336bce6e8SMatthew G. Knepley   PetscErrorCode ierr;
99436bce6e8SMatthew G. Knepley 
99536bce6e8SMatthew G. Knepley   PetscFunctionBegin;
9965cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
9974302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent, 2);
998c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
9994302210dSVaclav Hapla   PetscValidLogicalCollectiveEnum(viewer,datatype,4);
1000b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
10014302210dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr);
10024302210dSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr);
10034302210dSVaclav Hapla   ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);CHKERRQ(ierr);
100436bce6e8SMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
10057e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
10067e97c476SMatthew G. Knepley     size_t len;
10077e97c476SMatthew G. Knepley     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
1008729ab6d9SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
10097e97c476SMatthew G. Knepley   }
101036bce6e8SMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
1011729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
10124302210dSVaclav Hapla   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT));
1013080f144cSVaclav Hapla   if (has) {
1014080f144cSVaclav Hapla     PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
1015080f144cSVaclav Hapla   } else {
101660568592SMatthew G. Knepley     PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
1017080f144cSVaclav Hapla   }
1018729ab6d9SBarry Smith   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
1019729ab6d9SBarry Smith   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
1020729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
102160568592SMatthew G. Knepley   PetscStackCallHDF5(H5Oclose,(obj));
1022729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(dataspace));
10234302210dSVaclav Hapla   ierr = PetscFree(parentAbsPath);CHKERRQ(ierr);
102436bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
102536bce6e8SMatthew G. Knepley }
102636bce6e8SMatthew G. Knepley 
1027df863907SAlex Fikl /*@C
1028577e0e04SVaclav Hapla  PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name
1029577e0e04SVaclav Hapla 
1030343a20b2SVaclav Hapla   Collective
1031343a20b2SVaclav Hapla 
1032577e0e04SVaclav Hapla   Input Parameters:
1033577e0e04SVaclav Hapla + viewer   - The HDF5 viewer
1034577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1035577e0e04SVaclav Hapla . name     - The attribute name
1036577e0e04SVaclav Hapla . datatype - The attribute type
1037577e0e04SVaclav Hapla - value    - The attribute value
1038577e0e04SVaclav Hapla 
1039577e0e04SVaclav Hapla   Notes:
1040343a20b2SVaclav Hapla   This fails if the path current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1041577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
1042577e0e04SVaclav Hapla 
1043577e0e04SVaclav Hapla   Level: advanced
1044577e0e04SVaclav Hapla 
1045577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1046577e0e04SVaclav Hapla @*/
1047577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
1048577e0e04SVaclav Hapla {
1049577e0e04SVaclav Hapla   PetscErrorCode ierr;
1050577e0e04SVaclav Hapla 
1051577e0e04SVaclav Hapla   PetscFunctionBegin;
1052577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1053577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
1054577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
1055b7e81f0fSVaclav Hapla   PetscValidPointer(value,5);
1056577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
1057577e0e04SVaclav Hapla   ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
1058577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1059577e0e04SVaclav Hapla }
1060577e0e04SVaclav Hapla 
1061577e0e04SVaclav Hapla /*@C
1062b8ee85a1SVaclav Hapla  PetscViewerHDF5ReadAttribute - Read an attribute
1063ad0c61feSMatthew G. Knepley 
1064343a20b2SVaclav Hapla   Collective
1065343a20b2SVaclav Hapla 
1066ad0c61feSMatthew G. Knepley   Input Parameters:
1067ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer
10684302210dSVaclav Hapla . parent - The parent dataset/group name
1069ad0c61feSMatthew G. Knepley . name   - The attribute name
1070a2d6be1bSVaclav Hapla . datatype - The attribute type
1071a2d6be1bSVaclav Hapla - defaultValue - The pointer to the default value
1072ad0c61feSMatthew G. Knepley 
1073ad0c61feSMatthew G. Knepley   Output Parameter:
1074a2d6be1bSVaclav Hapla . value    - The pointer to the read HDF5 attribute value
1075ad0c61feSMatthew G. Knepley 
1076a2d6be1bSVaclav Hapla   Notes:
1077a2d6be1bSVaclav Hapla   If defaultValue is NULL and the attribute is not found, an error occurs.
1078a2d6be1bSVaclav Hapla   If defaultValue is not NULL and the attribute is not found, defaultValue is copied to value.
1079a2d6be1bSVaclav Hapla   The pointers defaultValue and value can be the same; for instance
1080a2d6be1bSVaclav Hapla $  flg = PETSC_FALSE;
1081a2d6be1bSVaclav Hapla $  ierr = PetscViewerHDF5ReadAttribute(viewer,name,"attr",PETSC_BOOL,&flg,&flg);CHKERRQ(ierr);
1082a2d6be1bSVaclav Hapla   is valid, but make sure the default value is initialized.
1083a2d6be1bSVaclav Hapla 
1084a2d6be1bSVaclav Hapla   If the datatype is PETSC_STRING, the output string is newly allocated so one must PetscFree() it when no longer needed.
1085708d4cb3SBarry Smith 
1086343a20b2SVaclav Hapla   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group. NULL means the current pushed group.
1087ad0c61feSMatthew G. Knepley 
1088343a20b2SVaclav Hapla   Level: advanced
10894302210dSVaclav Hapla 
1090577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1091ad0c61feSMatthew G. Knepley @*/
1092d70ec8faSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value)
1093ad0c61feSMatthew G. Knepley {
10944302210dSVaclav Hapla   char           *parentAbsPath;
1095a2d6be1bSVaclav Hapla   hid_t          h5, obj, attribute, dtype;
1096969748fdSVaclav Hapla   PetscBool      has;
1097ad0c61feSMatthew G. Knepley   PetscErrorCode ierr;
1098ad0c61feSMatthew G. Knepley 
1099ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
11005cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
11014302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent, 2);
1102c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
1103a2d6be1bSVaclav Hapla   if (defaultValue) PetscValidPointer(defaultValue, 5);
1104a2d6be1bSVaclav Hapla   PetscValidPointer(value, 6);
1105a2d6be1bSVaclav Hapla   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
11064302210dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr);
11074302210dSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL);CHKERRQ(ierr);
11084302210dSVaclav Hapla   if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);CHKERRQ(ierr);}
1109a2d6be1bSVaclav Hapla   if (!has) {
1110a2d6be1bSVaclav Hapla     if (defaultValue) {
1111a2d6be1bSVaclav Hapla       if (defaultValue != value) {
1112a2d6be1bSVaclav Hapla         if (datatype == PETSC_STRING) {
1113a2d6be1bSVaclav Hapla           ierr = PetscStrallocpy(*(char**)defaultValue, (char**)value);CHKERRQ(ierr);
1114a2d6be1bSVaclav Hapla         } else {
1115a2d6be1bSVaclav Hapla           size_t len;
1116a2d6be1bSVaclav Hapla           PetscStackCallHDF5Return(len,H5Tget_size,(dtype));
1117a2d6be1bSVaclav Hapla           ierr = PetscMemcpy(value, defaultValue, len);CHKERRQ(ierr);
1118a2d6be1bSVaclav Hapla         }
1119a2d6be1bSVaclav Hapla       }
1120a2d6be1bSVaclav Hapla       ierr = PetscFree(parentAbsPath);CHKERRQ(ierr);
1121a2d6be1bSVaclav Hapla       PetscFunctionReturn(0);
1122a2d6be1bSVaclav Hapla     } else SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name);
1123a2d6be1bSVaclav Hapla   }
1124ad0c61feSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
11254302210dSVaclav Hapla   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT));
112660568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
1127f0b58479SMatthew G. Knepley   if (datatype == PETSC_STRING) {
1128f0b58479SMatthew G. Knepley     size_t len;
1129a2d6be1bSVaclav Hapla     hid_t  atype;
1130e90c5075SMarek Pecha     PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
1131a2d6be1bSVaclav Hapla     PetscStackCallHDF5Return(len,H5Tget_size,(atype));
1132708d4cb3SBarry Smith     ierr = PetscMalloc((len+1) * sizeof(char), value);CHKERRQ(ierr);
1133708d4cb3SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
1134708d4cb3SBarry Smith     PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value));
1135708d4cb3SBarry Smith   } else {
113670efba33SBarry Smith     PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
1137708d4cb3SBarry Smith   }
1138729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
1139e90c5075SMarek Pecha   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
1140e90c5075SMarek Pecha   PetscStackCallHDF5(H5Oclose,(obj));
11414302210dSVaclav Hapla   ierr = PetscFree(parentAbsPath);CHKERRQ(ierr);
1142ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
1143ad0c61feSMatthew G. Knepley }
1144ad0c61feSMatthew G. Knepley 
1145577e0e04SVaclav Hapla /*@C
1146577e0e04SVaclav Hapla  PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name
1147577e0e04SVaclav Hapla 
1148343a20b2SVaclav Hapla   Collective
1149343a20b2SVaclav Hapla 
1150577e0e04SVaclav Hapla   Input Parameters:
1151577e0e04SVaclav Hapla + viewer   - The HDF5 viewer
1152577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1153577e0e04SVaclav Hapla . name     - The attribute name
1154577e0e04SVaclav Hapla - datatype - The attribute type
1155577e0e04SVaclav Hapla 
1156577e0e04SVaclav Hapla   Output Parameter:
1157577e0e04SVaclav Hapla . value    - The attribute value
1158577e0e04SVaclav Hapla 
1159577e0e04SVaclav Hapla   Notes:
1160577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1161577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
1162577e0e04SVaclav Hapla 
1163577e0e04SVaclav Hapla   Level: advanced
1164577e0e04SVaclav Hapla 
1165577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1166577e0e04SVaclav Hapla @*/
1167a2d6be1bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value)
1168577e0e04SVaclav Hapla {
1169577e0e04SVaclav Hapla   PetscErrorCode ierr;
1170577e0e04SVaclav Hapla 
1171577e0e04SVaclav Hapla   PetscFunctionBegin;
1172577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1173577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
1174577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
1175064a246eSJacob Faibussowitsch   PetscValidPointer(value, 6);
1176577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
1177a2d6be1bSVaclav Hapla   ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value);CHKERRQ(ierr);
1178577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1179577e0e04SVaclav Hapla }
1180577e0e04SVaclav Hapla 
1181a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
1182a75c3fd4SVaclav Hapla {
1183a75c3fd4SVaclav Hapla   htri_t exists;
1184a75c3fd4SVaclav Hapla   hid_t group;
1185a75c3fd4SVaclav Hapla 
1186a75c3fd4SVaclav Hapla   PetscFunctionBegin;
1187a75c3fd4SVaclav Hapla   PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT));
1188a75c3fd4SVaclav Hapla   if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT));
1189a75c3fd4SVaclav Hapla   if (!exists && createGroup) {
1190a75c3fd4SVaclav Hapla     PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1191a75c3fd4SVaclav Hapla     PetscStackCallHDF5(H5Gclose,(group));
1192a75c3fd4SVaclav Hapla     exists = PETSC_TRUE;
1193a75c3fd4SVaclav Hapla   }
1194a75c3fd4SVaclav Hapla   *exists_ = (PetscBool) exists;
1195a75c3fd4SVaclav Hapla   PetscFunctionReturn(0);
1196a75c3fd4SVaclav Hapla }
1197a75c3fd4SVaclav Hapla 
1198bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
11995cdeb108SMatthew G. Knepley {
120090e03537SVaclav Hapla   const char     rootGroupName[] = "/";
12015cdeb108SMatthew G. Knepley   hid_t          h5;
1202e5bf9ebcSVaclav Hapla   PetscBool      exists=PETSC_FALSE;
120315b861d2SVaclav Hapla   PetscInt       i;
120415b861d2SVaclav Hapla   int            n;
120585688ae2SVaclav Hapla   char           **hierarchy;
120685688ae2SVaclav Hapla   char           buf[PETSC_MAX_PATH_LEN]="";
12075cdeb108SMatthew G. Knepley   PetscErrorCode ierr;
12085cdeb108SMatthew G. Knepley 
12095cdeb108SMatthew G. Knepley   PetscFunctionBegin;
12105cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
121190e03537SVaclav Hapla   if (name) PetscValidCharPointer(name, 2);
121290e03537SVaclav Hapla   else name = rootGroupName;
1213ccd66a83SVaclav Hapla   if (has) {
1214064a246eSJacob Faibussowitsch     PetscValidBoolPointer(has, 4);
12155cdeb108SMatthew G. Knepley     *has = PETSC_FALSE;
1216ccd66a83SVaclav Hapla   }
1217ccd66a83SVaclav Hapla   if (otype) {
1218064a246eSJacob Faibussowitsch     PetscValidIntPointer(otype, 5);
121956cc0592SVaclav Hapla     *otype = H5O_TYPE_UNKNOWN;
1220ccd66a83SVaclav Hapla   }
1221c0cd0301SJed Brown   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
122285688ae2SVaclav Hapla 
122385688ae2SVaclav Hapla   /*
122485688ae2SVaclav Hapla      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
122585688ae2SVaclav Hapla      Hence, each of them needs to be tested separately:
122685688ae2SVaclav Hapla      1) whether it's a valid link
122785688ae2SVaclav Hapla      2) whether this link resolves to an object
122885688ae2SVaclav Hapla      See H5Oexists_by_name() documentation.
122985688ae2SVaclav Hapla   */
123085688ae2SVaclav Hapla   ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr);
123185688ae2SVaclav Hapla   if (!n) {
123285688ae2SVaclav Hapla     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1233ccd66a83SVaclav Hapla     if (has)   *has   = PETSC_TRUE;
1234ccd66a83SVaclav Hapla     if (otype) *otype = H5O_TYPE_GROUP;
123585688ae2SVaclav Hapla     ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
123685688ae2SVaclav Hapla     PetscFunctionReturn(0);
123785688ae2SVaclav Hapla   }
123885688ae2SVaclav Hapla   for (i=0; i<n; i++) {
123985688ae2SVaclav Hapla     ierr = PetscStrcat(buf,"/");CHKERRQ(ierr);
124085688ae2SVaclav Hapla     ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr);
1241a75c3fd4SVaclav Hapla     ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr);
1242a75c3fd4SVaclav Hapla     if (!exists) break;
124385688ae2SVaclav Hapla   }
124485688ae2SVaclav Hapla   ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
124585688ae2SVaclav Hapla 
124685688ae2SVaclav Hapla   /* If the object exists, get its type */
1247ccd66a83SVaclav Hapla   if (exists && otype) {
12485cdeb108SMatthew G. Knepley     H5O_info_t info;
12495cdeb108SMatthew G. Knepley 
125076276748SVaclav Hapla     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
125104633f7fSVaclav Hapla     PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT));
125256cc0592SVaclav Hapla     *otype = info.type;
12535cdeb108SMatthew G. Knepley   }
1254ccd66a83SVaclav Hapla   if (has) *has = exists;
12555cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
12565cdeb108SMatthew G. Knepley }
12575cdeb108SMatthew G. Knepley 
12584302210dSVaclav Hapla /*@C
12590a9f38efSVaclav Hapla  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
12600a9f38efSVaclav Hapla 
1261343a20b2SVaclav Hapla   Collective
1262343a20b2SVaclav Hapla 
12630a9f38efSVaclav Hapla   Input Parameters:
12644302210dSVaclav Hapla + viewer - The HDF5 viewer
12654302210dSVaclav Hapla - path - The group path
12660a9f38efSVaclav Hapla 
12670a9f38efSVaclav Hapla   Output Parameter:
12680a9f38efSVaclav Hapla . has - Flag for group existence
12690a9f38efSVaclav Hapla 
12700a9f38efSVaclav Hapla   Level: advanced
12710a9f38efSVaclav Hapla 
12724302210dSVaclav Hapla   Notes:
1273785c443eSVaclav Hapla   If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group.
1274785c443eSVaclav Hapla   So NULL or empty path means the current pushed group.
1275343a20b2SVaclav Hapla   If path exists but is not a group, PETSC_FALSE is returned.
12764302210dSVaclav Hapla 
127789e0ef10SVaclav Hapla .seealso: PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasDataset(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup(), PetscViewerHDF5OpenGroup()
12780a9f38efSVaclav Hapla @*/
12794302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has)
12800a9f38efSVaclav Hapla {
12810a9f38efSVaclav Hapla   H5O_type_t     type;
12824302210dSVaclav Hapla   char           *abspath;
12830a9f38efSVaclav Hapla   PetscErrorCode ierr;
12840a9f38efSVaclav Hapla 
12850a9f38efSVaclav Hapla   PetscFunctionBegin;
12860a9f38efSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
12874302210dSVaclav Hapla   if (path) PetscValidCharPointer(path,2);
12884302210dSVaclav Hapla   PetscValidBoolPointer(has,3);
12894302210dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);CHKERRQ(ierr);
12904302210dSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);CHKERRQ(ierr);
12914302210dSVaclav Hapla   *has = (PetscBool)(type == H5O_TYPE_GROUP);
12924302210dSVaclav Hapla   ierr = PetscFree(abspath);CHKERRQ(ierr);
12930a9f38efSVaclav Hapla   PetscFunctionReturn(0);
12940a9f38efSVaclav Hapla }
12950a9f38efSVaclav Hapla 
129689e0ef10SVaclav Hapla /*@C
129789e0ef10SVaclav Hapla  PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file
129889e0ef10SVaclav Hapla 
1299343a20b2SVaclav Hapla   Collective
1300343a20b2SVaclav Hapla 
130189e0ef10SVaclav Hapla   Input Parameters:
130289e0ef10SVaclav Hapla + viewer - The HDF5 viewer
130389e0ef10SVaclav Hapla - path - The dataset path
130489e0ef10SVaclav Hapla 
130589e0ef10SVaclav Hapla   Output Parameter:
130689e0ef10SVaclav Hapla . has - Flag whether dataset exists
130789e0ef10SVaclav Hapla 
130889e0ef10SVaclav Hapla   Level: advanced
130989e0ef10SVaclav Hapla 
131089e0ef10SVaclav Hapla   Notes:
1311343a20b2SVaclav Hapla   If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group.
131289e0ef10SVaclav Hapla   If path is NULL or empty, has is set to PETSC_FALSE.
1313343a20b2SVaclav Hapla   If path exists but is not a dataset, has is set to PETSC_FALSE as well.
131489e0ef10SVaclav Hapla 
131589e0ef10SVaclav Hapla .seealso: PetscViewerHDF5HasObject(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasGroup(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup()
131689e0ef10SVaclav Hapla @*/
131789e0ef10SVaclav Hapla PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has)
131889e0ef10SVaclav Hapla {
131989e0ef10SVaclav Hapla   H5O_type_t     type;
132089e0ef10SVaclav Hapla   char           *abspath;
132189e0ef10SVaclav Hapla   PetscErrorCode ierr;
132289e0ef10SVaclav Hapla 
132389e0ef10SVaclav Hapla   PetscFunctionBegin;
132489e0ef10SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
132589e0ef10SVaclav Hapla   if (path) PetscValidCharPointer(path,2);
132689e0ef10SVaclav Hapla   PetscValidBoolPointer(has,3);
132789e0ef10SVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);CHKERRQ(ierr);
132889e0ef10SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);CHKERRQ(ierr);
132989e0ef10SVaclav Hapla   *has = (PetscBool)(type == H5O_TYPE_DATASET);
133089e0ef10SVaclav Hapla   ierr = PetscFree(abspath);CHKERRQ(ierr);
133189e0ef10SVaclav Hapla   PetscFunctionReturn(0);
133289e0ef10SVaclav Hapla }
133389e0ef10SVaclav Hapla 
13340a9f38efSVaclav Hapla /*@
1335e3f143f7SVaclav Hapla  PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1336ecce1506SVaclav Hapla 
1337343a20b2SVaclav Hapla   Collective
1338343a20b2SVaclav Hapla 
1339ecce1506SVaclav Hapla   Input Parameters:
1340ecce1506SVaclav Hapla + viewer - The HDF5 viewer
1341ecce1506SVaclav Hapla - obj    - The named object
1342ecce1506SVaclav Hapla 
1343ecce1506SVaclav Hapla   Output Parameter:
134489e0ef10SVaclav Hapla . has    - Flag for dataset existence
1345ecce1506SVaclav Hapla 
1346e3f143f7SVaclav Hapla   Notes:
134789e0ef10SVaclav Hapla   If the object is unnamed, an error occurs.
1348343a20b2SVaclav Hapla   If the path current_group/object_name exists but is not a dataset, has is set to PETSC_FALSE as well.
1349e3f143f7SVaclav Hapla 
1350ecce1506SVaclav Hapla   Level: advanced
1351ecce1506SVaclav Hapla 
135289e0ef10SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasDataset(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup()
1353ecce1506SVaclav Hapla @*/
1354ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1355ecce1506SVaclav Hapla {
135689e0ef10SVaclav Hapla   size_t         len;
1357ecce1506SVaclav Hapla   PetscErrorCode ierr;
1358ecce1506SVaclav Hapla 
1359ecce1506SVaclav Hapla   PetscFunctionBegin;
1360c57a1a9aSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1361c57a1a9aSVaclav Hapla   PetscValidHeader(obj,2);
13624302210dSVaclav Hapla   PetscValidBoolPointer(has,3);
136389e0ef10SVaclav Hapla   ierr = PetscStrlen(obj->name, &len);CHKERRQ(ierr);
136489e0ef10SVaclav Hapla   if (!len) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
136589e0ef10SVaclav Hapla   ierr = PetscViewerHDF5HasDataset(viewer, obj->name, has);CHKERRQ(ierr);
1366ecce1506SVaclav Hapla   PetscFunctionReturn(0);
1367ecce1506SVaclav Hapla }
1368ecce1506SVaclav Hapla 
1369df863907SAlex Fikl /*@C
1370ece83b6cSVaclav Hapla  PetscViewerHDF5HasAttribute - Check whether an attribute exists
1371e7bdbf8eSMatthew G. Knepley 
1372343a20b2SVaclav Hapla   Collective
1373343a20b2SVaclav Hapla 
1374e7bdbf8eSMatthew G. Knepley   Input Parameters:
1375e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer
13764302210dSVaclav Hapla . parent - The parent dataset/group name
1377e7bdbf8eSMatthew G. Knepley - name   - The attribute name
1378e7bdbf8eSMatthew G. Knepley 
1379e7bdbf8eSMatthew G. Knepley   Output Parameter:
1380e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
1381e7bdbf8eSMatthew G. Knepley 
1382e7bdbf8eSMatthew G. Knepley   Level: advanced
1383e7bdbf8eSMatthew G. Knepley 
13844302210dSVaclav Hapla   Notes:
1385343a20b2SVaclav Hapla   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group. NULL means the current pushed group.
13864302210dSVaclav Hapla 
1387577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1388e7bdbf8eSMatthew G. Knepley @*/
13894302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1390e7bdbf8eSMatthew G. Knepley {
13914302210dSVaclav Hapla   char           *parentAbsPath;
1392e7bdbf8eSMatthew G. Knepley   PetscErrorCode ierr;
1393e7bdbf8eSMatthew G. Knepley 
1394e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
13955cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
13964302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent,2);
1397c57a1a9aSVaclav Hapla   PetscValidCharPointer(name,3);
13984302210dSVaclav Hapla   PetscValidBoolPointer(has,4);
13994302210dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr);
14004302210dSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL);CHKERRQ(ierr);
14014302210dSVaclav Hapla   if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has);CHKERRQ(ierr);}
14024302210dSVaclav Hapla   ierr = PetscFree(parentAbsPath);CHKERRQ(ierr);
140306db490cSVaclav Hapla   PetscFunctionReturn(0);
140406db490cSVaclav Hapla }
140506db490cSVaclav Hapla 
1406577e0e04SVaclav Hapla /*@C
1407577e0e04SVaclav Hapla  PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name
1408577e0e04SVaclav Hapla 
1409343a20b2SVaclav Hapla   Collective
1410343a20b2SVaclav Hapla 
1411577e0e04SVaclav Hapla   Input Parameters:
1412577e0e04SVaclav Hapla + viewer - The HDF5 viewer
1413577e0e04SVaclav Hapla . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1414577e0e04SVaclav Hapla - name   - The attribute name
1415577e0e04SVaclav Hapla 
1416577e0e04SVaclav Hapla   Output Parameter:
1417577e0e04SVaclav Hapla . has    - Flag for attribute existence
1418577e0e04SVaclav Hapla 
1419577e0e04SVaclav Hapla   Notes:
1420577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1421577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
1422577e0e04SVaclav Hapla 
1423577e0e04SVaclav Hapla   Level: advanced
1424577e0e04SVaclav Hapla 
1425577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1426577e0e04SVaclav Hapla @*/
1427577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1428577e0e04SVaclav Hapla {
1429577e0e04SVaclav Hapla   PetscErrorCode ierr;
1430577e0e04SVaclav Hapla 
1431577e0e04SVaclav Hapla   PetscFunctionBegin;
1432577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1433577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
1434577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
14354302210dSVaclav Hapla   PetscValidBoolPointer(has,4);
1436577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
1437577e0e04SVaclav Hapla   ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr);
1438577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1439577e0e04SVaclav Hapla }
1440577e0e04SVaclav Hapla 
144106db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
144206db490cSVaclav Hapla {
144306db490cSVaclav Hapla   hid_t          h5;
144406db490cSVaclav Hapla   htri_t         hhas;
144506db490cSVaclav Hapla   PetscErrorCode ierr;
144606db490cSVaclav Hapla 
144706db490cSVaclav Hapla   PetscFunctionBegin;
144806db490cSVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
14492f936e54SVaclav Hapla   PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT));
14505cdeb108SMatthew G. Knepley   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1451e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
1452e7bdbf8eSMatthew G. Knepley }
1453e7bdbf8eSMatthew G. Knepley 
1454a75e6a4aSMatthew G. Knepley /*
1455a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1456a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
1457a75e6a4aSMatthew G. Knepley */
1458d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1459a75e6a4aSMatthew G. Knepley 
1460a75e6a4aSMatthew G. Knepley /*@C
1461a75e6a4aSMatthew G. Knepley   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
1462a75e6a4aSMatthew G. Knepley 
1463d083f849SBarry Smith   Collective
1464a75e6a4aSMatthew G. Knepley 
1465a75e6a4aSMatthew G. Knepley   Input Parameter:
1466a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer
1467a75e6a4aSMatthew G. Knepley 
1468a75e6a4aSMatthew G. Knepley   Level: intermediate
1469a75e6a4aSMatthew G. Knepley 
1470a75e6a4aSMatthew G. Knepley   Options Database Keys:
1471a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name>
1472a75e6a4aSMatthew G. Knepley 
1473a75e6a4aSMatthew G. Knepley   Environmental variables:
1474a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME
1475a75e6a4aSMatthew G. Knepley 
1476a75e6a4aSMatthew G. Knepley   Notes:
1477a75e6a4aSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
1478a75e6a4aSMatthew G. Knepley   an error code.  The HDF5 PetscViewer is usually used in the form
1479a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1480a75e6a4aSMatthew G. Knepley 
1481a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1482a75e6a4aSMatthew G. Knepley @*/
1483a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1484a75e6a4aSMatthew G. Knepley {
1485a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
1486a75e6a4aSMatthew G. Knepley   PetscBool      flg;
1487a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
1488a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
1489a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
1490a75e6a4aSMatthew G. Knepley 
1491a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
1492908793a3SLisandro Dalcin   ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1493a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1494908793a3SLisandro Dalcin     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,NULL);
1495908793a3SLisandro Dalcin     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1496a75e6a4aSMatthew G. Knepley   }
149747435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1498908793a3SLisandro Dalcin   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1499a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
1500a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
15012cb5e1ccSBarry Smith     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1502a75e6a4aSMatthew G. Knepley     if (!flg) {
1503a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname,"output.h5");
15042cb5e1ccSBarry Smith       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1505a75e6a4aSMatthew G. Knepley     }
1506a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
15072cb5e1ccSBarry Smith     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1508a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
15092cb5e1ccSBarry Smith     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
151047435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
1511908793a3SLisandro Dalcin     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1512a75e6a4aSMatthew G. Knepley   }
1513a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
15142cb5e1ccSBarry Smith   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1515a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
1516a75e6a4aSMatthew G. Knepley }
1517