xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 534a8f05a7a8aff70dd8cfd53d9cd834400a8dbf)
1bc307625SVaclav Hapla #include <petsc/private/viewerimpl.h>
27d964842SVaclav Hapla #include <petsc/private/viewerhdf5impl.h>
3d70abbfaSBarry Smith #include <petscviewerhdf5.h>    /*I   "petscviewerhdf5.h"   I*/
45c6c1daeSBarry Smith 
5bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool*, H5O_type_t*);
606db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool*);
706db490cSVaclav Hapla 
86c132bc1SVaclav Hapla static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char objname[], char **fullpath)
96c132bc1SVaclav Hapla {
106c132bc1SVaclav Hapla   const char *group;
116c132bc1SVaclav Hapla   char buf[PETSC_MAX_PATH_LEN]="";
126c132bc1SVaclav Hapla   PetscErrorCode ierr;
136c132bc1SVaclav Hapla 
146c132bc1SVaclav Hapla   PetscFunctionBegin;
156c132bc1SVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr);
166c132bc1SVaclav Hapla   ierr = PetscStrcat(buf, group);CHKERRQ(ierr);
176c132bc1SVaclav Hapla   ierr = PetscStrcat(buf, "/");CHKERRQ(ierr);
186c132bc1SVaclav Hapla   ierr = PetscStrcat(buf, objname);CHKERRQ(ierr);
196c132bc1SVaclav Hapla   ierr = PetscStrallocpy(buf, fullpath);CHKERRQ(ierr);
206c132bc1SVaclav Hapla   PetscFunctionReturn(0);
216c132bc1SVaclav Hapla }
226c132bc1SVaclav Hapla 
23577e0e04SVaclav Hapla static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj)
24577e0e04SVaclav Hapla {
25577e0e04SVaclav Hapla   PetscBool has;
26577e0e04SVaclav Hapla   const char *group;
27577e0e04SVaclav Hapla   PetscErrorCode ierr;
28577e0e04SVaclav Hapla 
29577e0e04SVaclav Hapla   PetscFunctionBegin;
30577e0e04SVaclav Hapla   if (!obj->name) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
31577e0e04SVaclav Hapla   ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr);
32577e0e04SVaclav Hapla   if (!has) {
33577e0e04SVaclav Hapla     ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr);
34577e0e04SVaclav Hapla     SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) %s not stored in group %s", obj->name, group);
35577e0e04SVaclav Hapla   }
36577e0e04SVaclav Hapla   PetscFunctionReturn(0);
37577e0e04SVaclav Hapla }
38577e0e04SVaclav Hapla 
394416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v)
4082ea9b62SBarry Smith {
4182ea9b62SBarry Smith   PetscErrorCode   ierr;
42ee8b9732SVaclav Hapla   PetscBool        flg = PETSC_FALSE, set;
4382ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;
4482ea9b62SBarry Smith 
4582ea9b62SBarry Smith   PetscFunctionBegin;
4682ea9b62SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr);
4782ea9b62SBarry Smith   ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr);
489a0502c6SHåkon Strandenes   ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr);
49ee8b9732SVaclav Hapla   ierr = PetscOptionsBool("-viewer_hdf5_collective","Enable collective transfer mode","PetscViewerHDF5SetCollective",flg,&flg,&set);CHKERRQ(ierr);
50ee8b9732SVaclav Hapla   if (set) {ierr = PetscViewerHDF5SetCollective(v,flg);CHKERRQ(ierr);}
5182ea9b62SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
5282ea9b62SBarry Smith   PetscFunctionReturn(0);
5382ea9b62SBarry Smith }
5482ea9b62SBarry Smith 
551b793a25SVaclav Hapla static PetscErrorCode PetscViewerView_HDF5(PetscViewer v,PetscViewer viewer)
561b793a25SVaclav Hapla {
571b793a25SVaclav Hapla   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*)v->data;
5803fa8834SVaclav Hapla   PetscBool         flg;
591b793a25SVaclav Hapla   PetscErrorCode    ierr;
601b793a25SVaclav Hapla 
611b793a25SVaclav Hapla   PetscFunctionBegin;
621b793a25SVaclav Hapla   if (hdf5->filename) {
631b793a25SVaclav Hapla     ierr = PetscViewerASCIIPrintf(viewer,"Filename: %s\n",hdf5->filename);CHKERRQ(ierr);
641b793a25SVaclav Hapla   }
651b793a25SVaclav Hapla   ierr = PetscViewerASCIIPrintf(viewer,"Vectors with blocksize 1 saved as 2D datasets: %s\n",PetscBools[hdf5->basedimension2]);CHKERRQ(ierr);
661b793a25SVaclav Hapla   ierr = PetscViewerASCIIPrintf(viewer,"Enforce single precision storage: %s\n",PetscBools[hdf5->spoutput]);CHKERRQ(ierr);
6703fa8834SVaclav Hapla   ierr = PetscViewerHDF5GetCollective(v,&flg);CHKERRQ(ierr);
6803fa8834SVaclav Hapla   ierr = PetscViewerASCIIPrintf(viewer,"MPI-IO transfer mode: %s\n",flg ? "collective" : "independent");CHKERRQ(ierr);
691b793a25SVaclav Hapla   PetscFunctionReturn(0);
701b793a25SVaclav Hapla }
711b793a25SVaclav Hapla 
725c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
735c6c1daeSBarry Smith {
745c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
755c6c1daeSBarry Smith   PetscErrorCode   ierr;
765c6c1daeSBarry Smith 
775c6c1daeSBarry Smith   PetscFunctionBegin;
785c6c1daeSBarry Smith   ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);
79729ab6d9SBarry Smith   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
805c6c1daeSBarry Smith   PetscFunctionReturn(0);
815c6c1daeSBarry Smith }
825c6c1daeSBarry Smith 
839b1403beSVaclav Hapla static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
845c6c1daeSBarry Smith {
855c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
865c6c1daeSBarry Smith   PetscErrorCode   ierr;
875c6c1daeSBarry Smith 
885c6c1daeSBarry Smith   PetscFunctionBegin;
899c5072fbSVaclav Hapla   PetscStackCallHDF5(H5Pclose,(hdf5->dxpl_id));
905c6c1daeSBarry Smith   ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr);
915c6c1daeSBarry Smith   while (hdf5->groups) {
927d964842SVaclav Hapla     PetscViewerHDF5GroupList *tmp = hdf5->groups->next;
935c6c1daeSBarry Smith 
945c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups->name);CHKERRQ(ierr);
955c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups);CHKERRQ(ierr);
965c6c1daeSBarry Smith     hdf5->groups = tmp;
975c6c1daeSBarry Smith   }
985c6c1daeSBarry Smith   ierr = PetscFree(hdf5);CHKERRQ(ierr);
990b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
100d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr);
1010b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
102058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr);
103058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr);
1045c6c1daeSBarry Smith   PetscFunctionReturn(0);
1055c6c1daeSBarry Smith }
1065c6c1daeSBarry Smith 
1079b1403beSVaclav Hapla static PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
1085c6c1daeSBarry Smith {
1095c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1105c6c1daeSBarry Smith 
1115c6c1daeSBarry Smith   PetscFunctionBegin;
1125c6c1daeSBarry Smith   hdf5->btype = type;
1135c6c1daeSBarry Smith   PetscFunctionReturn(0);
1145c6c1daeSBarry Smith }
1155c6c1daeSBarry Smith 
1160b62783dSVaclav Hapla static PetscErrorCode  PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
1170b62783dSVaclav Hapla {
1180b62783dSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1190b62783dSVaclav Hapla 
1200b62783dSVaclav Hapla   PetscFunctionBegin;
1210b62783dSVaclav Hapla   *type = hdf5->btype;
1220b62783dSVaclav Hapla   PetscFunctionReturn(0);
1230b62783dSVaclav Hapla }
1240b62783dSVaclav Hapla 
1259b1403beSVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
12682ea9b62SBarry Smith {
12782ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
12882ea9b62SBarry Smith 
12982ea9b62SBarry Smith   PetscFunctionBegin;
13082ea9b62SBarry Smith   hdf5->basedimension2 = flg;
13182ea9b62SBarry Smith   PetscFunctionReturn(0);
13282ea9b62SBarry Smith }
13382ea9b62SBarry Smith 
134df863907SAlex Fikl /*@
13582ea9b62SBarry Smith      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
13682ea9b62SBarry Smith        dimension of 2.
13782ea9b62SBarry Smith 
13882ea9b62SBarry Smith     Logically Collective on PetscViewer
13982ea9b62SBarry Smith 
14082ea9b62SBarry Smith   Input Parameters:
14182ea9b62SBarry Smith +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
14282ea9b62SBarry 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
14382ea9b62SBarry Smith 
14482ea9b62SBarry Smith   Options Database:
14582ea9b62SBarry 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
14682ea9b62SBarry Smith 
14782ea9b62SBarry Smith 
14895452b02SPatrick Sanan   Notes:
14995452b02SPatrick 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
15082ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
15182ea9b62SBarry Smith 
15282ea9b62SBarry Smith   Level: intermediate
15382ea9b62SBarry Smith 
15482ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
15582ea9b62SBarry Smith 
15682ea9b62SBarry Smith @*/
15782ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
15882ea9b62SBarry Smith {
15982ea9b62SBarry Smith   PetscErrorCode ierr;
16082ea9b62SBarry Smith 
16182ea9b62SBarry Smith   PetscFunctionBegin;
16282ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
16382ea9b62SBarry Smith   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
16482ea9b62SBarry Smith   PetscFunctionReturn(0);
16582ea9b62SBarry Smith }
16682ea9b62SBarry Smith 
167df863907SAlex Fikl /*@
16882ea9b62SBarry Smith      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
16982ea9b62SBarry Smith        dimension of 2.
17082ea9b62SBarry Smith 
17182ea9b62SBarry Smith     Logically Collective on PetscViewer
17282ea9b62SBarry Smith 
17382ea9b62SBarry Smith   Input Parameter:
17482ea9b62SBarry Smith .  viewer - the PetscViewer, must be of type HDF5
17582ea9b62SBarry Smith 
17682ea9b62SBarry Smith   Output Parameter:
17782ea9b62SBarry 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
17882ea9b62SBarry Smith 
17995452b02SPatrick Sanan   Notes:
18095452b02SPatrick 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
18182ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
18282ea9b62SBarry Smith 
18382ea9b62SBarry Smith   Level: intermediate
18482ea9b62SBarry Smith 
18582ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
18682ea9b62SBarry Smith 
18782ea9b62SBarry Smith @*/
18882ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
18982ea9b62SBarry Smith {
19082ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
19182ea9b62SBarry Smith 
19282ea9b62SBarry Smith   PetscFunctionBegin;
19382ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
19482ea9b62SBarry Smith   *flg = hdf5->basedimension2;
19582ea9b62SBarry Smith   PetscFunctionReturn(0);
19682ea9b62SBarry Smith }
19782ea9b62SBarry Smith 
1989b1403beSVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
1999a0502c6SHåkon Strandenes {
2009a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2019a0502c6SHåkon Strandenes 
2029a0502c6SHåkon Strandenes   PetscFunctionBegin;
2039a0502c6SHåkon Strandenes   hdf5->spoutput = flg;
2049a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2059a0502c6SHåkon Strandenes }
2069a0502c6SHåkon Strandenes 
207df863907SAlex Fikl /*@
2089a0502c6SHåkon Strandenes      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
2099a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2109a0502c6SHåkon Strandenes 
2119a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2129a0502c6SHåkon Strandenes 
2139a0502c6SHåkon Strandenes   Input Parameters:
2149a0502c6SHåkon Strandenes +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
2159a0502c6SHåkon Strandenes -  flg - if PETSC_TRUE the data will be written to disk with single precision
2169a0502c6SHåkon Strandenes 
2179a0502c6SHåkon Strandenes   Options Database:
2189a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
2199a0502c6SHåkon Strandenes 
2209a0502c6SHåkon Strandenes 
22195452b02SPatrick Sanan   Notes:
22295452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2239a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2249a0502c6SHåkon Strandenes 
2259a0502c6SHåkon Strandenes   Level: intermediate
2269a0502c6SHåkon Strandenes 
2279a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2289a0502c6SHåkon Strandenes           PetscReal
2299a0502c6SHåkon Strandenes 
2309a0502c6SHåkon Strandenes @*/
2319a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
2329a0502c6SHåkon Strandenes {
2339a0502c6SHåkon Strandenes   PetscErrorCode ierr;
2349a0502c6SHåkon Strandenes 
2359a0502c6SHåkon Strandenes   PetscFunctionBegin;
2369a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2379a0502c6SHåkon Strandenes   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
2389a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2399a0502c6SHåkon Strandenes }
2409a0502c6SHåkon Strandenes 
241df863907SAlex Fikl /*@
2429a0502c6SHåkon Strandenes      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
2439a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2449a0502c6SHåkon Strandenes 
2459a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2469a0502c6SHåkon Strandenes 
2479a0502c6SHåkon Strandenes   Input Parameter:
2489a0502c6SHåkon Strandenes .  viewer - the PetscViewer, must be of type HDF5
2499a0502c6SHåkon Strandenes 
2509a0502c6SHåkon Strandenes   Output Parameter:
2519a0502c6SHåkon Strandenes .  flg - if PETSC_TRUE the data will be written to disk with single precision
2529a0502c6SHåkon Strandenes 
25395452b02SPatrick Sanan   Notes:
25495452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2559a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2569a0502c6SHåkon Strandenes 
2579a0502c6SHåkon Strandenes   Level: intermediate
2589a0502c6SHåkon Strandenes 
2599a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2609a0502c6SHåkon Strandenes           PetscReal
2619a0502c6SHåkon Strandenes 
2629a0502c6SHåkon Strandenes @*/
2639a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
2649a0502c6SHåkon Strandenes {
2659a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2669a0502c6SHåkon Strandenes 
2679a0502c6SHåkon Strandenes   PetscFunctionBegin;
2689a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2699a0502c6SHåkon Strandenes   *flg = hdf5->spoutput;
2709a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2719a0502c6SHåkon Strandenes }
2729a0502c6SHåkon Strandenes 
273ee8b9732SVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg)
274ee8b9732SVaclav Hapla {
275ee8b9732SVaclav Hapla   PetscFunctionBegin;
276ee8b9732SVaclav Hapla   /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions
277ee8b9732SVaclav Hapla      - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */
278ee8b9732SVaclav Hapla #if H5_VERSION_GE(1,10,3)
279ee8b9732SVaclav Hapla   {
280ee8b9732SVaclav Hapla     PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
281ee8b9732SVaclav Hapla     PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT));
282ee8b9732SVaclav Hapla   }
283ee8b9732SVaclav Hapla #else
284ee8b9732SVaclav Hapla   if (flg) {
285ee8b9732SVaclav Hapla     PetscErrorCode ierr;
286ee8b9732SVaclav Hapla     ierr = PetscPrintf(PetscObjectComm((PetscObject)viewer), "Warning: PetscViewerHDF5SetCollective(viewer,PETSC_TRUE) is ignored for HDF5 versions prior to 1.10.3\n");CHKERRQ(ierr);
287ee8b9732SVaclav Hapla   }
288ee8b9732SVaclav Hapla #endif
289ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
290ee8b9732SVaclav Hapla }
291ee8b9732SVaclav Hapla 
292ee8b9732SVaclav Hapla /*@
293ee8b9732SVaclav Hapla   PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes.
294ee8b9732SVaclav Hapla 
295ee8b9732SVaclav Hapla   Logically Collective; flg must contain common value
296ee8b9732SVaclav Hapla 
297ee8b9732SVaclav Hapla   Input Parameters:
298ee8b9732SVaclav Hapla + viewer - the PetscViewer; if it is not hdf5 then this command is ignored
299ee8b9732SVaclav Hapla - flg - PETSC_TRUE for collective mode; PETSC_FALSE for independent mode (default)
300ee8b9732SVaclav Hapla 
301ee8b9732SVaclav Hapla   Options Database:
302ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers
303ee8b9732SVaclav Hapla 
304ee8b9732SVaclav Hapla   Notes:
305ee8b9732SVaclav Hapla   Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better.
306ee8b9732SVaclav Hapla   However, this works correctly only since HDF5 1.10.3; hence, we ignore this setting for older versions.
307ee8b9732SVaclav Hapla 
308ee8b9732SVaclav Hapla   Developer notes:
309ee8b9732SVaclav Hapla   In the HDF5 layer, PETSC_TRUE / PETSC_FALSE means H5Pset_dxpl_mpio() is called with H5FD_MPIO_COLLECTIVE / H5FD_MPIO_INDEPENDENT, respectively.
310ee8b9732SVaclav Hapla   This in turn means use of MPI_File_{read,write}_all /  MPI_File_{read,write} in the MPI-IO layer, respectively.
311ee8b9732SVaclav Hapla   See HDF5 documentation and MPI-IO documentation for details.
312ee8b9732SVaclav Hapla 
313ee8b9732SVaclav Hapla   Level: intermediate
314ee8b9732SVaclav Hapla 
315ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5GetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open()
316ee8b9732SVaclav Hapla 
317ee8b9732SVaclav Hapla @*/
318ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg)
319ee8b9732SVaclav Hapla {
320ee8b9732SVaclav Hapla   PetscErrorCode ierr;
321ee8b9732SVaclav Hapla 
322ee8b9732SVaclav Hapla   PetscFunctionBegin;
323ee8b9732SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
324ee8b9732SVaclav Hapla   PetscValidLogicalCollectiveBool(viewer,flg,2);
325ee8b9732SVaclav Hapla   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
326ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
327ee8b9732SVaclav Hapla }
328ee8b9732SVaclav Hapla 
329ee8b9732SVaclav Hapla static PetscErrorCode  PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg)
330ee8b9732SVaclav Hapla {
331ee8b9732SVaclav Hapla   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*) viewer->data;
332ee8b9732SVaclav Hapla   H5FD_mpio_xfer_t  mode;
333ee8b9732SVaclav Hapla 
334ee8b9732SVaclav Hapla   PetscFunctionBegin;
335ee8b9732SVaclav Hapla   PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode));
336ee8b9732SVaclav Hapla   *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE;
337ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
338ee8b9732SVaclav Hapla }
339ee8b9732SVaclav Hapla 
340ee8b9732SVaclav Hapla /*@
341ee8b9732SVaclav Hapla   PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes.
342ee8b9732SVaclav Hapla 
343ee8b9732SVaclav Hapla   Not Collective
344ee8b9732SVaclav Hapla 
345ee8b9732SVaclav Hapla   Input Parameters:
346ee8b9732SVaclav Hapla . viewer - the HDF5 PetscViewer
347ee8b9732SVaclav Hapla 
348ee8b9732SVaclav Hapla   Output Parameters:
349ee8b9732SVaclav Hapla . flg - the flag
350ee8b9732SVaclav Hapla 
351ee8b9732SVaclav Hapla   Level: intermediate
352ee8b9732SVaclav Hapla 
353ee8b9732SVaclav Hapla   Notes:
354ee8b9732SVaclav Hapla   This setting works correctly only since HDF5 1.10.3. For older versions, PETSC_FALSE will be always returned.
355ee8b9732SVaclav Hapla   For more details, see PetscViewerHDF5SetCollective().
356ee8b9732SVaclav Hapla 
357ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open()
358ee8b9732SVaclav Hapla 
359ee8b9732SVaclav Hapla @*/
360ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg)
361ee8b9732SVaclav Hapla {
362ee8b9732SVaclav Hapla   PetscErrorCode ierr;
363ee8b9732SVaclav Hapla 
364ee8b9732SVaclav Hapla   PetscFunctionBegin;
365ee8b9732SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
366*534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg,2);
367ee8b9732SVaclav Hapla 
368ee8b9732SVaclav Hapla   ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));CHKERRQ(ierr);
369ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
370ee8b9732SVaclav Hapla }
371ee8b9732SVaclav Hapla 
3729b1403beSVaclav Hapla static PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
3735c6c1daeSBarry Smith {
3745c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3755c6c1daeSBarry Smith   hid_t             plist_id;
3765c6c1daeSBarry Smith   PetscErrorCode    ierr;
3775c6c1daeSBarry Smith 
3785c6c1daeSBarry Smith   PetscFunctionBegin;
379f683cc58SToby Isaac   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
380f683cc58SToby Isaac   if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);}
3815c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
3825c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
383729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
384d28bfa55SVaclav Hapla   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL));
3855c6c1daeSBarry Smith   /* Create or open the file collectively */
3865c6c1daeSBarry Smith   switch (hdf5->btype) {
3875c6c1daeSBarry Smith   case FILE_MODE_READ:
388729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
3895c6c1daeSBarry Smith     break;
3905c6c1daeSBarry Smith   case FILE_MODE_APPEND:
391729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
3925c6c1daeSBarry Smith     break;
3935c6c1daeSBarry Smith   case FILE_MODE_WRITE:
394729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
3955c6c1daeSBarry Smith     break;
3965c6c1daeSBarry Smith   default:
3975c6c1daeSBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
3985c6c1daeSBarry Smith   }
3995c6c1daeSBarry Smith   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
400729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
4015c6c1daeSBarry Smith   PetscFunctionReturn(0);
4025c6c1daeSBarry Smith }
4035c6c1daeSBarry Smith 
404d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
405d1232d7fSToby Isaac {
406d1232d7fSToby Isaac   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
407d1232d7fSToby Isaac 
408d1232d7fSToby Isaac   PetscFunctionBegin;
409d1232d7fSToby Isaac   *name = vhdf5->filename;
410d1232d7fSToby Isaac   PetscFunctionReturn(0);
411d1232d7fSToby Isaac }
412d1232d7fSToby Isaac 
413b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
414b723ab35SVaclav Hapla {
4155dc64a97SVaclav Hapla   /*
416b723ab35SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
417b723ab35SVaclav Hapla   PetscErrorCode   ierr;
4185dc64a97SVaclav Hapla   */
419b723ab35SVaclav Hapla 
420b723ab35SVaclav Hapla   PetscFunctionBegin;
421b723ab35SVaclav Hapla   PetscFunctionReturn(0);
422b723ab35SVaclav Hapla }
423b723ab35SVaclav Hapla 
4248556b5ebSBarry Smith /*MC
4258556b5ebSBarry Smith    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
4268556b5ebSBarry Smith 
4278556b5ebSBarry Smith 
4288556b5ebSBarry Smith .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
4298556b5ebSBarry Smith            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
4308556b5ebSBarry Smith            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
4318556b5ebSBarry Smith            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
4328556b5ebSBarry Smith 
4331b266c99SBarry Smith   Level: beginner
4348556b5ebSBarry Smith M*/
435d1232d7fSToby Isaac 
4368cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
4375c6c1daeSBarry Smith {
4385c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
4395c6c1daeSBarry Smith   PetscErrorCode   ierr;
4405c6c1daeSBarry Smith 
4415c6c1daeSBarry Smith   PetscFunctionBegin;
442b00a9115SJed Brown   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
4435c6c1daeSBarry Smith 
4445c6c1daeSBarry Smith   v->data                = (void*) hdf5;
4455c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
44682ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
447b723ab35SVaclav Hapla   v->ops->setup          = PetscViewerSetUp_HDF5;
4481b793a25SVaclav Hapla   v->ops->view           = PetscViewerView_HDF5;
4495c6c1daeSBarry Smith   v->ops->flush          = 0;
4505c6c1daeSBarry Smith   hdf5->btype            = (PetscFileMode) -1;
4515c6c1daeSBarry Smith   hdf5->filename         = 0;
4525c6c1daeSBarry Smith   hdf5->timestep         = -1;
4530298fd71SBarry Smith   hdf5->groups           = NULL;
4545c6c1daeSBarry Smith 
4559c5072fbSVaclav Hapla   PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER));
4569c5072fbSVaclav Hapla 
4570b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
458d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr);
4590b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
4600b62783dSVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr);
46182ea9b62SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
4629a0502c6SHåkon Strandenes   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr);
463ee8b9732SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);CHKERRQ(ierr);
464ee8b9732SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);CHKERRQ(ierr);
4655c6c1daeSBarry Smith   PetscFunctionReturn(0);
4665c6c1daeSBarry Smith }
4675c6c1daeSBarry Smith 
4685c6c1daeSBarry Smith /*@C
4695c6c1daeSBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
4705c6c1daeSBarry Smith 
471d083f849SBarry Smith    Collective
4725c6c1daeSBarry Smith 
4735c6c1daeSBarry Smith    Input Parameters:
4745c6c1daeSBarry Smith +  comm - MPI communicator
4755c6c1daeSBarry Smith .  name - name of file
4765c6c1daeSBarry Smith -  type - type of file
4775c6c1daeSBarry Smith $    FILE_MODE_WRITE - create new file for binary output
4785c6c1daeSBarry Smith $    FILE_MODE_READ - open existing file for binary input
4795c6c1daeSBarry Smith $    FILE_MODE_APPEND - open existing file for binary output
4805c6c1daeSBarry Smith 
4815c6c1daeSBarry Smith    Output Parameter:
4825c6c1daeSBarry Smith .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
4835c6c1daeSBarry Smith 
48482ea9b62SBarry Smith   Options Database:
485a2b725a8SWilliam 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
486a2b725a8SWilliam Gropp -  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
48782ea9b62SBarry Smith 
4885c6c1daeSBarry Smith    Level: beginner
4895c6c1daeSBarry Smith 
4905c6c1daeSBarry Smith    Note:
4915c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
4925c6c1daeSBarry Smith 
4935c6c1daeSBarry Smith 
4946a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
4959a0502c6SHåkon Strandenes           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
496a56f64adSBarry Smith           MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
4975c6c1daeSBarry Smith @*/
4985c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
4995c6c1daeSBarry Smith {
5005c6c1daeSBarry Smith   PetscErrorCode ierr;
5015c6c1daeSBarry Smith 
5025c6c1daeSBarry Smith   PetscFunctionBegin;
5035c6c1daeSBarry Smith   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
5045c6c1daeSBarry Smith   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
5055c6c1daeSBarry Smith   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
5065c6c1daeSBarry Smith   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
507b663d2d2SVaclav Hapla   ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr);
5085c6c1daeSBarry Smith   PetscFunctionReturn(0);
5095c6c1daeSBarry Smith }
5105c6c1daeSBarry Smith 
5115c6c1daeSBarry Smith /*@C
5125c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
5135c6c1daeSBarry Smith 
5145c6c1daeSBarry Smith   Not collective
5155c6c1daeSBarry Smith 
5165c6c1daeSBarry Smith   Input Parameter:
5175c6c1daeSBarry Smith . viewer - the PetscViewer
5185c6c1daeSBarry Smith 
5195c6c1daeSBarry Smith   Output Parameter:
5205c6c1daeSBarry Smith . file_id - The file id
5215c6c1daeSBarry Smith 
5225c6c1daeSBarry Smith   Level: intermediate
5235c6c1daeSBarry Smith 
5245c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open()
5255c6c1daeSBarry Smith @*/
5265c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
5275c6c1daeSBarry Smith {
5285c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5295c6c1daeSBarry Smith 
5305c6c1daeSBarry Smith   PetscFunctionBegin;
5315c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5325c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
5335c6c1daeSBarry Smith   PetscFunctionReturn(0);
5345c6c1daeSBarry Smith }
5355c6c1daeSBarry Smith 
5365c6c1daeSBarry Smith /*@C
5375c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
5385c6c1daeSBarry Smith 
5395c6c1daeSBarry Smith   Not collective
5405c6c1daeSBarry Smith 
5415c6c1daeSBarry Smith   Input Parameters:
5425c6c1daeSBarry Smith + viewer - the PetscViewer
5435c6c1daeSBarry Smith - name - The group name
5445c6c1daeSBarry Smith 
5455c6c1daeSBarry Smith   Level: intermediate
5465c6c1daeSBarry Smith 
547874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
5485c6c1daeSBarry Smith @*/
5495c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
5505c6c1daeSBarry Smith {
5515c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5527d964842SVaclav Hapla   PetscViewerHDF5GroupList *groupNode;
5535c6c1daeSBarry Smith   PetscErrorCode   ierr;
5545c6c1daeSBarry Smith 
5555c6c1daeSBarry Smith   PetscFunctionBegin;
5565c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5575c6c1daeSBarry Smith   PetscValidCharPointer(name,2);
55895dccacaSBarry Smith   ierr = PetscNew(&groupNode);CHKERRQ(ierr);
5595c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
560a297a907SKarl Rupp 
5615c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
5625c6c1daeSBarry Smith   hdf5->groups    = groupNode;
5635c6c1daeSBarry Smith   PetscFunctionReturn(0);
5645c6c1daeSBarry Smith }
5655c6c1daeSBarry Smith 
5663ef9c667SSatish Balay /*@
5675c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
5685c6c1daeSBarry Smith 
5695c6c1daeSBarry Smith   Not collective
5705c6c1daeSBarry Smith 
5715c6c1daeSBarry Smith   Input Parameter:
5725c6c1daeSBarry Smith . viewer - the PetscViewer
5735c6c1daeSBarry Smith 
5745c6c1daeSBarry Smith   Level: intermediate
5755c6c1daeSBarry Smith 
576874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
5775c6c1daeSBarry Smith @*/
5785c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
5795c6c1daeSBarry Smith {
5805c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5817d964842SVaclav Hapla   PetscViewerHDF5GroupList *groupNode;
5825c6c1daeSBarry Smith   PetscErrorCode   ierr;
5835c6c1daeSBarry Smith 
5845c6c1daeSBarry Smith   PetscFunctionBegin;
5855c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
58682f516ccSBarry Smith   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
5875c6c1daeSBarry Smith   groupNode    = hdf5->groups;
5885c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
5895c6c1daeSBarry Smith   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
5905c6c1daeSBarry Smith   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
5915c6c1daeSBarry Smith   PetscFunctionReturn(0);
5925c6c1daeSBarry Smith }
5935c6c1daeSBarry Smith 
5945c6c1daeSBarry Smith /*@C
595874270d9SVaclav Hapla   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup().
596874270d9SVaclav Hapla   If none has been assigned, returns NULL.
5975c6c1daeSBarry Smith 
5985c6c1daeSBarry Smith   Not collective
5995c6c1daeSBarry Smith 
6005c6c1daeSBarry Smith   Input Parameter:
6015c6c1daeSBarry Smith . viewer - the PetscViewer
6025c6c1daeSBarry Smith 
6035c6c1daeSBarry Smith   Output Parameter:
6045c6c1daeSBarry Smith . name - The group name
6055c6c1daeSBarry Smith 
6065c6c1daeSBarry Smith   Level: intermediate
6075c6c1daeSBarry Smith 
608874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
6095c6c1daeSBarry Smith @*/
6105c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
6115c6c1daeSBarry Smith {
6125c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
6135c6c1daeSBarry Smith 
6145c6c1daeSBarry Smith   PetscFunctionBegin;
6155c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
616c959eef4SJed Brown   PetscValidPointer(name,2);
617a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
6180298fd71SBarry Smith   else *name = NULL;
6195c6c1daeSBarry Smith   PetscFunctionReturn(0);
6205c6c1daeSBarry Smith }
6215c6c1daeSBarry Smith 
6228c557b5aSVaclav Hapla /*@
623874270d9SVaclav Hapla   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(),
624874270d9SVaclav Hapla   and return this group's ID and file ID.
625874270d9SVaclav Hapla   If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID.
626874270d9SVaclav Hapla 
627874270d9SVaclav Hapla   Not collective
628874270d9SVaclav Hapla 
629874270d9SVaclav Hapla   Input Parameter:
630874270d9SVaclav Hapla . viewer - the PetscViewer
631874270d9SVaclav Hapla 
632874270d9SVaclav Hapla   Output Parameter:
633874270d9SVaclav Hapla + fileId - The HDF5 file ID
634874270d9SVaclav Hapla - groupId - The HDF5 group ID
635874270d9SVaclav Hapla 
636e74616beSVaclav Hapla   Notes:
637e74616beSVaclav Hapla   If the viewer is writable, the group is created if it doesn't exist yet.
638e74616beSVaclav Hapla 
639874270d9SVaclav Hapla   Level: intermediate
640874270d9SVaclav Hapla 
641874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
642874270d9SVaclav Hapla @*/
64354dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
64454dbf706SVaclav Hapla {
64590e03537SVaclav Hapla   hid_t          file_id;
64690e03537SVaclav Hapla   H5O_type_t     type;
64754dbf706SVaclav Hapla   const char     *groupName = NULL;
648e74616beSVaclav Hapla   PetscBool      create;
64954dbf706SVaclav Hapla   PetscErrorCode ierr;
65054dbf706SVaclav Hapla 
65154dbf706SVaclav Hapla   PetscFunctionBegin;
652e74616beSVaclav Hapla   ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr);
65354dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr);
65454dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr);
655e74616beSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr);
65690e03537SVaclav 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);
65790e03537SVaclav Hapla   PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT));
65854dbf706SVaclav Hapla   *fileId  = file_id;
65954dbf706SVaclav Hapla   PetscFunctionReturn(0);
66054dbf706SVaclav Hapla }
66154dbf706SVaclav Hapla 
6623ef9c667SSatish Balay /*@
6635c6c1daeSBarry Smith   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
6645c6c1daeSBarry Smith 
6655c6c1daeSBarry Smith   Not collective
6665c6c1daeSBarry Smith 
6675c6c1daeSBarry Smith   Input Parameter:
6685c6c1daeSBarry Smith . viewer - the PetscViewer
6695c6c1daeSBarry Smith 
6705c6c1daeSBarry Smith   Level: intermediate
6715c6c1daeSBarry Smith 
6725c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
6735c6c1daeSBarry Smith @*/
6745c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
6755c6c1daeSBarry Smith {
6765c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
6775c6c1daeSBarry Smith 
6785c6c1daeSBarry Smith   PetscFunctionBegin;
6795c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
6805c6c1daeSBarry Smith   ++hdf5->timestep;
6815c6c1daeSBarry Smith   PetscFunctionReturn(0);
6825c6c1daeSBarry Smith }
6835c6c1daeSBarry Smith 
6843ef9c667SSatish Balay /*@
6855c6c1daeSBarry Smith   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
6865c6c1daeSBarry Smith   of -1 disables blocking with timesteps.
6875c6c1daeSBarry Smith 
6885c6c1daeSBarry Smith   Not collective
6895c6c1daeSBarry Smith 
6905c6c1daeSBarry Smith   Input Parameters:
6915c6c1daeSBarry Smith + viewer - the PetscViewer
6925c6c1daeSBarry Smith - timestep - The timestep number
6935c6c1daeSBarry Smith 
6945c6c1daeSBarry Smith   Level: intermediate
6955c6c1daeSBarry Smith 
6965c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
6975c6c1daeSBarry Smith @*/
6985c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
6995c6c1daeSBarry Smith {
7005c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7015c6c1daeSBarry Smith 
7025c6c1daeSBarry Smith   PetscFunctionBegin;
7035c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7045c6c1daeSBarry Smith   hdf5->timestep = timestep;
7055c6c1daeSBarry Smith   PetscFunctionReturn(0);
7065c6c1daeSBarry Smith }
7075c6c1daeSBarry Smith 
7083ef9c667SSatish Balay /*@
7095c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
7105c6c1daeSBarry Smith 
7115c6c1daeSBarry Smith   Not collective
7125c6c1daeSBarry Smith 
7135c6c1daeSBarry Smith   Input Parameter:
7145c6c1daeSBarry Smith . viewer - the PetscViewer
7155c6c1daeSBarry Smith 
7165c6c1daeSBarry Smith   Output Parameter:
7175c6c1daeSBarry Smith . timestep - The timestep number
7185c6c1daeSBarry Smith 
7195c6c1daeSBarry Smith   Level: intermediate
7205c6c1daeSBarry Smith 
7215c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
7225c6c1daeSBarry Smith @*/
7235c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
7245c6c1daeSBarry Smith {
7255c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7265c6c1daeSBarry Smith 
7275c6c1daeSBarry Smith   PetscFunctionBegin;
7285c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7295c6c1daeSBarry Smith   PetscValidPointer(timestep,2);
7305c6c1daeSBarry Smith   *timestep = hdf5->timestep;
7315c6c1daeSBarry Smith   PetscFunctionReturn(0);
7325c6c1daeSBarry Smith }
7335c6c1daeSBarry Smith 
73436bce6e8SMatthew G. Knepley /*@C
73536bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
73636bce6e8SMatthew G. Knepley 
73736bce6e8SMatthew G. Knepley   Not collective
73836bce6e8SMatthew G. Knepley 
73936bce6e8SMatthew G. Knepley   Input Parameter:
74036bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
74136bce6e8SMatthew G. Knepley 
74236bce6e8SMatthew G. Knepley   Output Parameter:
74336bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
74436bce6e8SMatthew G. Knepley 
74536bce6e8SMatthew G. Knepley   Level: advanced
74636bce6e8SMatthew G. Knepley 
74736bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
74836bce6e8SMatthew G. Knepley @*/
74936bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
75036bce6e8SMatthew G. Knepley {
75136bce6e8SMatthew G. Knepley   PetscFunctionBegin;
75236bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
75336bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
75436bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_LLONG;
75536bce6e8SMatthew G. Knepley #else
75636bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_INT;
75736bce6e8SMatthew G. Knepley #endif
75836bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
75936bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
76036bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
761c9450970SBarry Smith   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_INT;
762de3ba262SVaclav Hapla   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_INT;
76336bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
76436bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
76536bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
7667e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
76736bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
76836bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
76936bce6e8SMatthew G. Knepley }
77036bce6e8SMatthew G. Knepley 
77136bce6e8SMatthew G. Knepley /*@C
77236bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
77336bce6e8SMatthew G. Knepley 
77436bce6e8SMatthew G. Knepley   Not collective
77536bce6e8SMatthew G. Knepley 
77636bce6e8SMatthew G. Knepley   Input Parameter:
77736bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
77836bce6e8SMatthew G. Knepley 
77936bce6e8SMatthew G. Knepley   Output Parameter:
78036bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
78136bce6e8SMatthew G. Knepley 
78236bce6e8SMatthew G. Knepley   Level: advanced
78336bce6e8SMatthew G. Knepley 
78436bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
78536bce6e8SMatthew G. Knepley @*/
78636bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
78736bce6e8SMatthew G. Knepley {
78836bce6e8SMatthew G. Knepley   PetscFunctionBegin;
78936bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
79036bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
79136bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
79236bce6e8SMatthew G. Knepley #else
79336bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
79436bce6e8SMatthew G. Knepley #endif
79536bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
79636bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
79736bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
79836bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
79936bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
80036bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
8017e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
80236bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
80336bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
80436bce6e8SMatthew G. Knepley }
80536bce6e8SMatthew G. Knepley 
806df863907SAlex Fikl /*@C
807b8ee85a1SVaclav Hapla  PetscViewerHDF5WriteAttribute - Write an attribute
80836bce6e8SMatthew G. Knepley 
80936bce6e8SMatthew G. Knepley   Input Parameters:
81036bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer
81157d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
81236bce6e8SMatthew G. Knepley . name   - The attribute name
81336bce6e8SMatthew G. Knepley . datatype - The attribute type
81436bce6e8SMatthew G. Knepley - value    - The attribute value
81536bce6e8SMatthew G. Knepley 
81636bce6e8SMatthew G. Knepley   Level: advanced
81736bce6e8SMatthew G. Knepley 
818577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
81936bce6e8SMatthew G. Knepley @*/
82057d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value)
82136bce6e8SMatthew G. Knepley {
82257d22f7dSVaclav Hapla   char           *parent;
82360568592SMatthew G. Knepley   hid_t          h5, dataspace, obj, attribute, dtype;
824080f144cSVaclav Hapla   PetscBool      has;
82536bce6e8SMatthew G. Knepley   PetscErrorCode ierr;
82636bce6e8SMatthew G. Knepley 
82736bce6e8SMatthew G. Knepley   PetscFunctionBegin;
8285cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
82957d22f7dSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset, 2);
830c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
831b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
83257d22f7dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
833bb286ee1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr);
834b67695ecSVaclav Hapla   ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);
83536bce6e8SMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
8367e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
8377e97c476SMatthew G. Knepley     size_t len;
8387e97c476SMatthew G. Knepley     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
839729ab6d9SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
8407e97c476SMatthew G. Knepley   }
84136bce6e8SMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
842729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
84360568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
844080f144cSVaclav Hapla   if (has) {
845080f144cSVaclav Hapla     PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
846080f144cSVaclav Hapla   } else {
84760568592SMatthew G. Knepley     PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
848080f144cSVaclav Hapla   }
849729ab6d9SBarry Smith   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
850729ab6d9SBarry Smith   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
851729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
85260568592SMatthew G. Knepley   PetscStackCallHDF5(H5Oclose,(obj));
853729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(dataspace));
85457d22f7dSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
85536bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
85636bce6e8SMatthew G. Knepley }
85736bce6e8SMatthew G. Knepley 
858df863907SAlex Fikl /*@C
859577e0e04SVaclav Hapla  PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name
860577e0e04SVaclav Hapla 
861577e0e04SVaclav Hapla   Input Parameters:
862577e0e04SVaclav Hapla + viewer   - The HDF5 viewer
863577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
864577e0e04SVaclav Hapla . name     - The attribute name
865577e0e04SVaclav Hapla . datatype - The attribute type
866577e0e04SVaclav Hapla - value    - The attribute value
867577e0e04SVaclav Hapla 
868577e0e04SVaclav Hapla   Notes:
869577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
870577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
871577e0e04SVaclav Hapla 
872577e0e04SVaclav Hapla   Level: advanced
873577e0e04SVaclav Hapla 
874577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
875577e0e04SVaclav Hapla @*/
876577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
877577e0e04SVaclav Hapla {
878577e0e04SVaclav Hapla   PetscErrorCode ierr;
879577e0e04SVaclav Hapla 
880577e0e04SVaclav Hapla   PetscFunctionBegin;
881577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
882577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
883577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
884b7e81f0fSVaclav Hapla   PetscValidPointer(value,5);
885577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
886577e0e04SVaclav Hapla   ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
887577e0e04SVaclav Hapla   PetscFunctionReturn(0);
888577e0e04SVaclav Hapla }
889577e0e04SVaclav Hapla 
890577e0e04SVaclav Hapla /*@C
891b8ee85a1SVaclav Hapla  PetscViewerHDF5ReadAttribute - Read an attribute
892ad0c61feSMatthew G. Knepley 
893ad0c61feSMatthew G. Knepley   Input Parameters:
894ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer
89557d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
896ad0c61feSMatthew G. Knepley . name   - The attribute name
897ad0c61feSMatthew G. Knepley - datatype - The attribute type
898ad0c61feSMatthew G. Knepley 
899ad0c61feSMatthew G. Knepley   Output Parameter:
900ad0c61feSMatthew G. Knepley . value    - The attribute value
901ad0c61feSMatthew G. Knepley 
902ad0c61feSMatthew G. Knepley   Level: advanced
903ad0c61feSMatthew G. Knepley 
904577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
905ad0c61feSMatthew G. Knepley @*/
90657d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value)
907ad0c61feSMatthew G. Knepley {
90857d22f7dSVaclav Hapla   char           *parent;
909f0b58479SMatthew G. Knepley   hid_t          h5, obj, attribute, atype, dtype;
910969748fdSVaclav Hapla   PetscBool      has;
911ad0c61feSMatthew G. Knepley   PetscErrorCode ierr;
912ad0c61feSMatthew G. Knepley 
913ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
9145cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
91557d22f7dSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset, 2);
916c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
917b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
91857d22f7dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
919969748fdSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr);
920969748fdSVaclav Hapla   if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);}
921969748fdSVaclav Hapla   if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name);
922ad0c61feSMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
923ad0c61feSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
92460568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
92560568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
926f0b58479SMatthew G. Knepley   if (datatype == PETSC_STRING) {
927f0b58479SMatthew G. Knepley     size_t len;
928e90c5075SMarek Pecha     PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
92915b861d2SVaclav Hapla     PetscStackCall("H5Tget_size",len = H5Tget_size(atype));
930f0b58479SMatthew G. Knepley     PetscStackCallHDF5(H5Tclose,(atype));
931f0b58479SMatthew G. Knepley     ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr);
932f0b58479SMatthew G. Knepley   }
93370efba33SBarry Smith   PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
934729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
935e90c5075SMarek Pecha   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
936e90c5075SMarek Pecha   PetscStackCallHDF5(H5Oclose,(obj));
93757d22f7dSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
938ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
939ad0c61feSMatthew G. Knepley }
940ad0c61feSMatthew G. Knepley 
941577e0e04SVaclav Hapla /*@C
942577e0e04SVaclav Hapla  PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name
943577e0e04SVaclav Hapla 
944577e0e04SVaclav Hapla   Input Parameters:
945577e0e04SVaclav Hapla + viewer   - The HDF5 viewer
946577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
947577e0e04SVaclav Hapla . name     - The attribute name
948577e0e04SVaclav Hapla - datatype - The attribute type
949577e0e04SVaclav Hapla 
950577e0e04SVaclav Hapla   Output Parameter:
951577e0e04SVaclav Hapla . value    - The attribute value
952577e0e04SVaclav Hapla 
953577e0e04SVaclav Hapla   Notes:
954577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
955577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
956577e0e04SVaclav Hapla 
957577e0e04SVaclav Hapla   Level: advanced
958577e0e04SVaclav Hapla 
959577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
960577e0e04SVaclav Hapla @*/
961577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value)
962577e0e04SVaclav Hapla {
963577e0e04SVaclav Hapla   PetscErrorCode ierr;
964577e0e04SVaclav Hapla 
965577e0e04SVaclav Hapla   PetscFunctionBegin;
966577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
967577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
968577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
969b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
970577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
971577e0e04SVaclav Hapla   ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
972577e0e04SVaclav Hapla   PetscFunctionReturn(0);
973577e0e04SVaclav Hapla }
974577e0e04SVaclav Hapla 
975a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
976a75c3fd4SVaclav Hapla {
977a75c3fd4SVaclav Hapla   htri_t exists;
978a75c3fd4SVaclav Hapla   hid_t group;
979a75c3fd4SVaclav Hapla 
980a75c3fd4SVaclav Hapla   PetscFunctionBegin;
981a75c3fd4SVaclav Hapla   PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT));
982a75c3fd4SVaclav Hapla   if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT));
983a75c3fd4SVaclav Hapla   if (!exists && createGroup) {
984a75c3fd4SVaclav Hapla     PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
985a75c3fd4SVaclav Hapla     PetscStackCallHDF5(H5Gclose,(group));
986a75c3fd4SVaclav Hapla     exists = PETSC_TRUE;
987a75c3fd4SVaclav Hapla   }
988a75c3fd4SVaclav Hapla   *exists_ = (PetscBool) exists;
989a75c3fd4SVaclav Hapla   PetscFunctionReturn(0);
990a75c3fd4SVaclav Hapla }
991a75c3fd4SVaclav Hapla 
992bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
9935cdeb108SMatthew G. Knepley {
99490e03537SVaclav Hapla   const char     rootGroupName[] = "/";
9955cdeb108SMatthew G. Knepley   hid_t          h5;
996e5bf9ebcSVaclav Hapla   PetscBool      exists=PETSC_FALSE;
99715b861d2SVaclav Hapla   PetscInt       i;
99815b861d2SVaclav Hapla   int            n;
99985688ae2SVaclav Hapla   char           **hierarchy;
100085688ae2SVaclav Hapla   char           buf[PETSC_MAX_PATH_LEN]="";
10015cdeb108SMatthew G. Knepley   PetscErrorCode ierr;
10025cdeb108SMatthew G. Knepley 
10035cdeb108SMatthew G. Knepley   PetscFunctionBegin;
10045cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
100590e03537SVaclav Hapla   if (name) PetscValidCharPointer(name, 2);
100690e03537SVaclav Hapla   else name = rootGroupName;
1007ccd66a83SVaclav Hapla   if (has) {
100856cc0592SVaclav Hapla     PetscValidIntPointer(has, 3);
10095cdeb108SMatthew G. Knepley     *has = PETSC_FALSE;
1010ccd66a83SVaclav Hapla   }
1011ccd66a83SVaclav Hapla   if (otype) {
1012ccd66a83SVaclav Hapla     PetscValidIntPointer(otype, 4);
101356cc0592SVaclav Hapla     *otype = H5O_TYPE_UNKNOWN;
1014ccd66a83SVaclav Hapla   }
1015c0cd0301SJed Brown   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
101685688ae2SVaclav Hapla 
101785688ae2SVaclav Hapla   /*
101885688ae2SVaclav Hapla      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
101985688ae2SVaclav Hapla      Hence, each of them needs to be tested separately:
102085688ae2SVaclav Hapla      1) whether it's a valid link
102185688ae2SVaclav Hapla      2) whether this link resolves to an object
102285688ae2SVaclav Hapla      See H5Oexists_by_name() documentation.
102385688ae2SVaclav Hapla   */
102485688ae2SVaclav Hapla   ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr);
102585688ae2SVaclav Hapla   if (!n) {
102685688ae2SVaclav Hapla     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1027ccd66a83SVaclav Hapla     if (has)   *has   = PETSC_TRUE;
1028ccd66a83SVaclav Hapla     if (otype) *otype = H5O_TYPE_GROUP;
102985688ae2SVaclav Hapla     ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
103085688ae2SVaclav Hapla     PetscFunctionReturn(0);
103185688ae2SVaclav Hapla   }
103285688ae2SVaclav Hapla   for (i=0; i<n; i++) {
103385688ae2SVaclav Hapla     ierr = PetscStrcat(buf,"/");CHKERRQ(ierr);
103485688ae2SVaclav Hapla     ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr);
1035a75c3fd4SVaclav Hapla     ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr);
1036a75c3fd4SVaclav Hapla     if (!exists) break;
103785688ae2SVaclav Hapla   }
103885688ae2SVaclav Hapla   ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
103985688ae2SVaclav Hapla 
104085688ae2SVaclav Hapla   /* If the object exists, get its type */
1041ccd66a83SVaclav Hapla   if (exists && otype) {
10425cdeb108SMatthew G. Knepley     H5O_info_t info;
10435cdeb108SMatthew G. Knepley 
104476276748SVaclav Hapla     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
104504633f7fSVaclav Hapla     PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT));
104656cc0592SVaclav Hapla     *otype = info.type;
10475cdeb108SMatthew G. Knepley   }
1048ccd66a83SVaclav Hapla   if (has) *has = exists;
10495cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
10505cdeb108SMatthew G. Knepley }
10515cdeb108SMatthew G. Knepley 
1052ecce1506SVaclav Hapla /*@
10530a9f38efSVaclav Hapla  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
10540a9f38efSVaclav Hapla 
10550a9f38efSVaclav Hapla   Input Parameters:
10560a9f38efSVaclav Hapla . viewer - The HDF5 viewer
10570a9f38efSVaclav Hapla 
10580a9f38efSVaclav Hapla   Output Parameter:
10590a9f38efSVaclav Hapla . has    - Flag for group existence
10600a9f38efSVaclav Hapla 
10610a9f38efSVaclav Hapla   Notes:
10620a9f38efSVaclav Hapla   If the path exists but is not a group, this returns PETSC_FALSE as well.
10630a9f38efSVaclav Hapla 
10640a9f38efSVaclav Hapla   Level: advanced
10650a9f38efSVaclav Hapla 
10660a9f38efSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup()
10670a9f38efSVaclav Hapla @*/
10680a9f38efSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has)
10690a9f38efSVaclav Hapla {
10700a9f38efSVaclav Hapla   H5O_type_t type;
10710a9f38efSVaclav Hapla   const char *name;
10720a9f38efSVaclav Hapla   PetscErrorCode ierr;
10730a9f38efSVaclav Hapla 
10740a9f38efSVaclav Hapla   PetscFunctionBegin;
10750a9f38efSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
10760a9f38efSVaclav Hapla   PetscValidIntPointer(has,2);
10770a9f38efSVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr);
10780a9f38efSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr);
10790a9f38efSVaclav Hapla   *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE;
10800a9f38efSVaclav Hapla   PetscFunctionReturn(0);
10810a9f38efSVaclav Hapla }
10820a9f38efSVaclav Hapla 
10830a9f38efSVaclav Hapla /*@
1084e3f143f7SVaclav Hapla  PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1085ecce1506SVaclav Hapla 
1086ecce1506SVaclav Hapla   Input Parameters:
1087ecce1506SVaclav Hapla + viewer - The HDF5 viewer
1088ecce1506SVaclav Hapla - obj    - The named object
1089ecce1506SVaclav Hapla 
1090ecce1506SVaclav Hapla   Output Parameter:
1091ecce1506SVaclav Hapla . has    - Flag for dataset existence; PETSC_FALSE for unnamed object
1092ecce1506SVaclav Hapla 
1093e3f143f7SVaclav Hapla   Notes:
1094e3f143f7SVaclav Hapla   If the path exists but is not a dataset, this returns PETSC_FALSE as well.
1095e3f143f7SVaclav Hapla 
1096ecce1506SVaclav Hapla   Level: advanced
1097ecce1506SVaclav Hapla 
1098e3f143f7SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1099ecce1506SVaclav Hapla @*/
1100ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1101ecce1506SVaclav Hapla {
110256cc0592SVaclav Hapla   H5O_type_t type;
11036c132bc1SVaclav Hapla   char *path;
1104ecce1506SVaclav Hapla   PetscErrorCode ierr;
1105ecce1506SVaclav Hapla 
1106ecce1506SVaclav Hapla   PetscFunctionBegin;
1107c57a1a9aSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1108c57a1a9aSVaclav Hapla   PetscValidHeader(obj,2);
1109c57a1a9aSVaclav Hapla   PetscValidIntPointer(has,3);
1110ecce1506SVaclav Hapla   *has = PETSC_FALSE;
1111e3f143f7SVaclav Hapla   if (!obj->name) PetscFunctionReturn(0);
11126c132bc1SVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr);
11136c132bc1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr);
111456cc0592SVaclav Hapla   *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE;
11156c132bc1SVaclav Hapla   ierr = PetscFree(path);CHKERRQ(ierr);
1116ecce1506SVaclav Hapla   PetscFunctionReturn(0);
1117ecce1506SVaclav Hapla }
1118ecce1506SVaclav Hapla 
1119df863907SAlex Fikl /*@C
1120ece83b6cSVaclav Hapla  PetscViewerHDF5HasAttribute - Check whether an attribute exists
1121e7bdbf8eSMatthew G. Knepley 
1122e7bdbf8eSMatthew G. Knepley   Input Parameters:
1123e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer
112410e69e8fSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
1125e7bdbf8eSMatthew G. Knepley - name   - The attribute name
1126e7bdbf8eSMatthew G. Knepley 
1127e7bdbf8eSMatthew G. Knepley   Output Parameter:
1128e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
1129e7bdbf8eSMatthew G. Knepley 
1130e7bdbf8eSMatthew G. Knepley   Level: advanced
1131e7bdbf8eSMatthew G. Knepley 
1132577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1133e7bdbf8eSMatthew G. Knepley @*/
113410e69e8fSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has)
1135e7bdbf8eSMatthew G. Knepley {
113610e69e8fSVaclav Hapla   char           *parent;
1137e7bdbf8eSMatthew G. Knepley   PetscErrorCode ierr;
1138e7bdbf8eSMatthew G. Knepley 
1139e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
11405cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
114110e69e8fSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset,2);
1142c57a1a9aSVaclav Hapla   PetscValidCharPointer(name,3);
1143c57a1a9aSVaclav Hapla   PetscValidIntPointer(has,4);
114410e69e8fSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
1145bb286ee1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr);
114610e69e8fSVaclav Hapla   if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);}
114710e69e8fSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
114806db490cSVaclav Hapla   PetscFunctionReturn(0);
114906db490cSVaclav Hapla }
115006db490cSVaclav Hapla 
1151577e0e04SVaclav Hapla /*@C
1152577e0e04SVaclav Hapla  PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name
1153577e0e04SVaclav Hapla 
1154577e0e04SVaclav Hapla   Input Parameters:
1155577e0e04SVaclav Hapla + viewer - The HDF5 viewer
1156577e0e04SVaclav Hapla . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1157577e0e04SVaclav Hapla - name   - The attribute name
1158577e0e04SVaclav Hapla 
1159577e0e04SVaclav Hapla   Output Parameter:
1160577e0e04SVaclav Hapla . has    - Flag for attribute existence
1161577e0e04SVaclav Hapla 
1162577e0e04SVaclav Hapla   Notes:
1163577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1164577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
1165577e0e04SVaclav Hapla 
1166577e0e04SVaclav Hapla   Level: advanced
1167577e0e04SVaclav Hapla 
1168577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1169577e0e04SVaclav Hapla @*/
1170577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1171577e0e04SVaclav Hapla {
1172577e0e04SVaclav Hapla   PetscErrorCode ierr;
1173577e0e04SVaclav Hapla 
1174577e0e04SVaclav Hapla   PetscFunctionBegin;
1175577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1176577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
1177577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
1178577e0e04SVaclav Hapla   PetscValidIntPointer(has,4);
1179577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
1180577e0e04SVaclav Hapla   ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr);
1181577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1182577e0e04SVaclav Hapla }
1183577e0e04SVaclav Hapla 
118406db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
118506db490cSVaclav Hapla {
118606db490cSVaclav Hapla   hid_t          h5;
118706db490cSVaclav Hapla   htri_t         hhas;
118806db490cSVaclav Hapla   PetscErrorCode ierr;
118906db490cSVaclav Hapla 
119006db490cSVaclav Hapla   PetscFunctionBegin;
119106db490cSVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
11922f936e54SVaclav Hapla   PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT));
11935cdeb108SMatthew G. Knepley   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1194e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
1195e7bdbf8eSMatthew G. Knepley }
1196e7bdbf8eSMatthew G. Knepley 
1197a75e6a4aSMatthew G. Knepley /*
1198a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1199a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
1200a75e6a4aSMatthew G. Knepley */
1201d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1202a75e6a4aSMatthew G. Knepley 
1203a75e6a4aSMatthew G. Knepley /*@C
1204a75e6a4aSMatthew G. Knepley   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
1205a75e6a4aSMatthew G. Knepley 
1206d083f849SBarry Smith   Collective
1207a75e6a4aSMatthew G. Knepley 
1208a75e6a4aSMatthew G. Knepley   Input Parameter:
1209a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer
1210a75e6a4aSMatthew G. Knepley 
1211a75e6a4aSMatthew G. Knepley   Level: intermediate
1212a75e6a4aSMatthew G. Knepley 
1213a75e6a4aSMatthew G. Knepley   Options Database Keys:
1214a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name>
1215a75e6a4aSMatthew G. Knepley 
1216a75e6a4aSMatthew G. Knepley   Environmental variables:
1217a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME
1218a75e6a4aSMatthew G. Knepley 
1219a75e6a4aSMatthew G. Knepley   Notes:
1220a75e6a4aSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
1221a75e6a4aSMatthew G. Knepley   an error code.  The HDF5 PetscViewer is usually used in the form
1222a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1223a75e6a4aSMatthew G. Knepley 
1224a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1225a75e6a4aSMatthew G. Knepley @*/
1226a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1227a75e6a4aSMatthew G. Knepley {
1228a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
1229a75e6a4aSMatthew G. Knepley   PetscBool      flg;
1230a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
1231a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
1232a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
1233a75e6a4aSMatthew G. Knepley 
1234a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
1235a75e6a4aSMatthew G. Knepley   ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1236a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
123712801b39SBarry Smith     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
1238a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1239a75e6a4aSMatthew G. Knepley   }
124047435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1241a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1242a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
1243a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
1244a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1245a75e6a4aSMatthew G. Knepley     if (!flg) {
1246a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname,"output.h5");
1247a75e6a4aSMatthew G. Knepley       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1248a75e6a4aSMatthew G. Knepley     }
1249a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
1250a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1251a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
1252a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
125347435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
1254a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1255a75e6a4aSMatthew G. Knepley   }
1256a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
1257a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1258a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
1259a75e6a4aSMatthew G. Knepley }
1260