xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 76d59af26514faa333c8049ee56c7c266f42dc45)
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 */
278c796909eSBarry Smith #if H5_VERSION_GE(1,10,3) && defined(H5_HAVE_PARALLEL)
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;
286c796909eSBarry 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);
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.
30653effed7SBarry 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.
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 {
331c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL)
332ee8b9732SVaclav Hapla   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*) viewer->data;
333ee8b9732SVaclav Hapla   H5FD_mpio_xfer_t  mode;
334c796909eSBarry Smith #endif
335ee8b9732SVaclav Hapla 
336ee8b9732SVaclav Hapla   PetscFunctionBegin;
337c796909eSBarry Smith #if !defined(H5_HAVE_PARALLEL)
338c796909eSBarry Smith   *flg = PETSC_FALSE;
339c796909eSBarry Smith #else
340ee8b9732SVaclav Hapla   PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode));
341ee8b9732SVaclav Hapla   *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE;
342c796909eSBarry Smith #endif
343ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
344ee8b9732SVaclav Hapla }
345ee8b9732SVaclav Hapla 
346ee8b9732SVaclav Hapla /*@
347ee8b9732SVaclav Hapla   PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes.
348ee8b9732SVaclav Hapla 
349ee8b9732SVaclav Hapla   Not Collective
350ee8b9732SVaclav Hapla 
351ee8b9732SVaclav Hapla   Input Parameters:
352ee8b9732SVaclav Hapla . viewer - the HDF5 PetscViewer
353ee8b9732SVaclav Hapla 
354ee8b9732SVaclav Hapla   Output Parameters:
355ee8b9732SVaclav Hapla . flg - the flag
356ee8b9732SVaclav Hapla 
357ee8b9732SVaclav Hapla   Level: intermediate
358ee8b9732SVaclav Hapla 
359ee8b9732SVaclav Hapla   Notes:
360c796909eSBarry 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.
361ee8b9732SVaclav Hapla   For more details, see PetscViewerHDF5SetCollective().
362ee8b9732SVaclav Hapla 
363ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open()
364ee8b9732SVaclav Hapla 
365ee8b9732SVaclav Hapla @*/
366ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg)
367ee8b9732SVaclav Hapla {
368ee8b9732SVaclav Hapla   PetscErrorCode ierr;
369ee8b9732SVaclav Hapla 
370ee8b9732SVaclav Hapla   PetscFunctionBegin;
371ee8b9732SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
372534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg,2);
373ee8b9732SVaclav Hapla 
374ee8b9732SVaclav Hapla   ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));CHKERRQ(ierr);
375ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
376ee8b9732SVaclav Hapla }
377ee8b9732SVaclav Hapla 
3789b1403beSVaclav Hapla static PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
3795c6c1daeSBarry Smith {
3805c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3815c6c1daeSBarry Smith   hid_t             plist_id;
3825c6c1daeSBarry Smith   PetscErrorCode    ierr;
3835c6c1daeSBarry Smith 
3845c6c1daeSBarry Smith   PetscFunctionBegin;
385f683cc58SToby Isaac   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
386f683cc58SToby Isaac   if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);}
3875c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
3885c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
389729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
390c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL)
391d28bfa55SVaclav Hapla   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL));
392c796909eSBarry Smith #endif
3935c6c1daeSBarry Smith   /* Create or open the file collectively */
3945c6c1daeSBarry Smith   switch (hdf5->btype) {
3955c6c1daeSBarry Smith   case FILE_MODE_READ:
396729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
3975c6c1daeSBarry Smith     break;
3985c6c1daeSBarry Smith   case FILE_MODE_APPEND:
3997e4fd573SVaclav Hapla   case FILE_MODE_UPDATE:
4007e4fd573SVaclav Hapla   {
4017e4fd573SVaclav Hapla     PetscBool flg;
4027e4fd573SVaclav Hapla     ierr = PetscTestFile(hdf5->filename, 'r', &flg);CHKERRQ(ierr);
4037e4fd573SVaclav Hapla     if (flg) PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
4047e4fd573SVaclav Hapla     else     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id));
4055c6c1daeSBarry Smith     break;
4067e4fd573SVaclav Hapla   }
4075c6c1daeSBarry Smith   case FILE_MODE_WRITE:
408729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
4095c6c1daeSBarry Smith     break;
4107e4fd573SVaclav Hapla   case FILE_MODE_UNDEFINED:
4117e4fd573SVaclav Hapla     SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
4125c6c1daeSBarry Smith   default:
4137e4fd573SVaclav Hapla     SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP, "Unsupported file mode %s",PetscFileModes[hdf5->btype]);
4145c6c1daeSBarry Smith   }
4155c6c1daeSBarry Smith   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
416729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
4175c6c1daeSBarry Smith   PetscFunctionReturn(0);
4185c6c1daeSBarry Smith }
4195c6c1daeSBarry Smith 
420d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
421d1232d7fSToby Isaac {
422d1232d7fSToby Isaac   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
423d1232d7fSToby Isaac 
424d1232d7fSToby Isaac   PetscFunctionBegin;
425d1232d7fSToby Isaac   *name = vhdf5->filename;
426d1232d7fSToby Isaac   PetscFunctionReturn(0);
427d1232d7fSToby Isaac }
428d1232d7fSToby Isaac 
429b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
430b723ab35SVaclav Hapla {
4315dc64a97SVaclav Hapla   /*
432b723ab35SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
433b723ab35SVaclav Hapla   PetscErrorCode   ierr;
4345dc64a97SVaclav Hapla   */
435b723ab35SVaclav Hapla 
436b723ab35SVaclav Hapla   PetscFunctionBegin;
437b723ab35SVaclav Hapla   PetscFunctionReturn(0);
438b723ab35SVaclav Hapla }
439b723ab35SVaclav Hapla 
4408556b5ebSBarry Smith /*MC
4418556b5ebSBarry Smith    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
4428556b5ebSBarry Smith 
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;
45899335c34SVaclav Hapla #if !defined(H5_HAVE_PARALLEL)
45999335c34SVaclav Hapla   {
46099335c34SVaclav Hapla     PetscMPIInt size;
46199335c34SVaclav Hapla     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)v), &size);CHKERRMPI(ierr);
46299335c34SVaclav Hapla     if (size > 1) SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "Cannot use parallel HDF5 viewer since the given HDF5 does not support parallel I/O (H5_HAVE_PARALLEL is unset)");
46399335c34SVaclav Hapla   }
46499335c34SVaclav Hapla #endif
46599335c34SVaclav Hapla 
466b00a9115SJed Brown   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
4675c6c1daeSBarry Smith 
4685c6c1daeSBarry Smith   v->data                = (void*) hdf5;
4695c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
47082ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
471b723ab35SVaclav Hapla   v->ops->setup          = PetscViewerSetUp_HDF5;
4721b793a25SVaclav Hapla   v->ops->view           = PetscViewerView_HDF5;
473908793a3SLisandro Dalcin   v->ops->flush          = NULL;
4747e4fd573SVaclav Hapla   hdf5->btype            = FILE_MODE_UNDEFINED;
475908793a3SLisandro Dalcin   hdf5->filename         = NULL;
4765c6c1daeSBarry Smith   hdf5->timestep         = -1;
4770298fd71SBarry Smith   hdf5->groups           = NULL;
4785c6c1daeSBarry Smith 
4799c5072fbSVaclav Hapla   PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER));
4809c5072fbSVaclav Hapla 
4810b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
482d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr);
4830b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
4840b62783dSVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr);
48582ea9b62SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
4869a0502c6SHåkon Strandenes   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr);
487ee8b9732SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);CHKERRQ(ierr);
488ee8b9732SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);CHKERRQ(ierr);
4895c6c1daeSBarry Smith   PetscFunctionReturn(0);
4905c6c1daeSBarry Smith }
4915c6c1daeSBarry Smith 
4925c6c1daeSBarry Smith /*@C
4935c6c1daeSBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
4945c6c1daeSBarry Smith 
495d083f849SBarry Smith    Collective
4965c6c1daeSBarry Smith 
4975c6c1daeSBarry Smith    Input Parameters:
4985c6c1daeSBarry Smith +  comm - MPI communicator
4995c6c1daeSBarry Smith .  name - name of file
5005c6c1daeSBarry Smith -  type - type of file
5015c6c1daeSBarry Smith 
5025c6c1daeSBarry Smith    Output Parameter:
5035c6c1daeSBarry Smith .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
5045c6c1daeSBarry Smith 
50582ea9b62SBarry Smith   Options Database:
506a2b725a8SWilliam 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
507a2b725a8SWilliam Gropp -  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
50882ea9b62SBarry Smith 
5095c6c1daeSBarry Smith    Level: beginner
5105c6c1daeSBarry Smith 
5117e4fd573SVaclav Hapla    Notes:
5127e4fd573SVaclav Hapla    Reading is always available, regardless of the mode. Available modes are
5137e4fd573SVaclav Hapla +  FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY]
5147e4fd573SVaclav Hapla .  FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC]
5157e4fd573SVaclav 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]
5167e4fd573SVaclav Hapla -  FILE_MODE_UPDATE - same as FILE_MODE_APPEND
5177e4fd573SVaclav Hapla 
5187e4fd573SVaclav 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.
5197e4fd573SVaclav Hapla 
5205c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
5215c6c1daeSBarry Smith 
5225c6c1daeSBarry Smith 
5236a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
5249a0502c6SHåkon Strandenes           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
525a56f64adSBarry Smith           MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
5265c6c1daeSBarry Smith @*/
5275c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
5285c6c1daeSBarry Smith {
5295c6c1daeSBarry Smith   PetscErrorCode ierr;
5305c6c1daeSBarry Smith 
5315c6c1daeSBarry Smith   PetscFunctionBegin;
5325c6c1daeSBarry Smith   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
5335c6c1daeSBarry Smith   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
5345c6c1daeSBarry Smith   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
5355c6c1daeSBarry Smith   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
536b663d2d2SVaclav Hapla   ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr);
5375c6c1daeSBarry Smith   PetscFunctionReturn(0);
5385c6c1daeSBarry Smith }
5395c6c1daeSBarry Smith 
5405c6c1daeSBarry Smith /*@C
5415c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
5425c6c1daeSBarry Smith 
5435c6c1daeSBarry Smith   Not collective
5445c6c1daeSBarry Smith 
5455c6c1daeSBarry Smith   Input Parameter:
5465c6c1daeSBarry Smith . viewer - the PetscViewer
5475c6c1daeSBarry Smith 
5485c6c1daeSBarry Smith   Output Parameter:
5495c6c1daeSBarry Smith . file_id - The file id
5505c6c1daeSBarry Smith 
5515c6c1daeSBarry Smith   Level: intermediate
5525c6c1daeSBarry Smith 
5535c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open()
5545c6c1daeSBarry Smith @*/
5555c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
5565c6c1daeSBarry Smith {
5575c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5585c6c1daeSBarry Smith 
5595c6c1daeSBarry Smith   PetscFunctionBegin;
5605c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5615c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
5625c6c1daeSBarry Smith   PetscFunctionReturn(0);
5635c6c1daeSBarry Smith }
5645c6c1daeSBarry Smith 
5655c6c1daeSBarry Smith /*@C
5665c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
5675c6c1daeSBarry Smith 
5685c6c1daeSBarry Smith   Not collective
5695c6c1daeSBarry Smith 
5705c6c1daeSBarry Smith   Input Parameters:
5715c6c1daeSBarry Smith + viewer - the PetscViewer
5725c6c1daeSBarry Smith - name - The group name
5735c6c1daeSBarry Smith 
5745c6c1daeSBarry Smith   Level: intermediate
5755c6c1daeSBarry Smith 
576820fc2d1SVaclav Hapla   Note: The group name being NULL, empty string, or a sequence of all slashes (e.g. "///") is always internally stored as NULL and interpreted as "/".
577820fc2d1SVaclav Hapla 
578874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
5795c6c1daeSBarry Smith @*/
580be7aa430SVaclav Hapla PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[])
5815c6c1daeSBarry Smith {
5825c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5837d964842SVaclav Hapla   PetscViewerHDF5GroupList *groupNode;
5845c6c1daeSBarry Smith   PetscErrorCode   ierr;
5855c6c1daeSBarry Smith 
5865c6c1daeSBarry Smith   PetscFunctionBegin;
5875c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
588820fc2d1SVaclav Hapla   if (name) PetscValidCharPointer(name,2);
589820fc2d1SVaclav Hapla   if (name && name[0]) {
590820fc2d1SVaclav Hapla      size_t i,len;
591820fc2d1SVaclav Hapla      ierr = PetscStrlen(name, &len);CHKERRQ(ierr);
592820fc2d1SVaclav Hapla      for (i=0; i<len; i++) if (name[i] != '/') break;
593820fc2d1SVaclav Hapla      if (i == len) name = NULL;
594820fc2d1SVaclav Hapla   } else name = NULL;
59595dccacaSBarry Smith   ierr = PetscNew(&groupNode);CHKERRQ(ierr);
5965c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
5975c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
5985c6c1daeSBarry Smith   hdf5->groups    = groupNode;
5995c6c1daeSBarry Smith   PetscFunctionReturn(0);
6005c6c1daeSBarry Smith }
6015c6c1daeSBarry Smith 
6023ef9c667SSatish Balay /*@
6035c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
6045c6c1daeSBarry Smith 
6055c6c1daeSBarry Smith   Not collective
6065c6c1daeSBarry Smith 
6075c6c1daeSBarry Smith   Input Parameter:
6085c6c1daeSBarry Smith . viewer - the PetscViewer
6095c6c1daeSBarry Smith 
6105c6c1daeSBarry Smith   Level: intermediate
6115c6c1daeSBarry Smith 
612874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
6135c6c1daeSBarry Smith @*/
6145c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
6155c6c1daeSBarry Smith {
6165c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
6177d964842SVaclav Hapla   PetscViewerHDF5GroupList *groupNode;
6185c6c1daeSBarry Smith   PetscErrorCode   ierr;
6195c6c1daeSBarry Smith 
6205c6c1daeSBarry Smith   PetscFunctionBegin;
6215c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
62282f516ccSBarry Smith   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
6235c6c1daeSBarry Smith   groupNode    = hdf5->groups;
6245c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
6255c6c1daeSBarry Smith   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
6265c6c1daeSBarry Smith   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
6275c6c1daeSBarry Smith   PetscFunctionReturn(0);
6285c6c1daeSBarry Smith }
6295c6c1daeSBarry Smith 
6305c6c1daeSBarry Smith /*@C
631874270d9SVaclav Hapla   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup().
632874270d9SVaclav Hapla   If none has been assigned, returns NULL.
6335c6c1daeSBarry Smith 
6345c6c1daeSBarry Smith   Not collective
6355c6c1daeSBarry Smith 
6365c6c1daeSBarry Smith   Input Parameter:
6375c6c1daeSBarry Smith . viewer - the PetscViewer
6385c6c1daeSBarry Smith 
6395c6c1daeSBarry Smith   Output Parameter:
6405c6c1daeSBarry Smith . name - The group name
6415c6c1daeSBarry Smith 
6425c6c1daeSBarry Smith   Level: intermediate
6435c6c1daeSBarry Smith 
644874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
6455c6c1daeSBarry Smith @*/
646be7aa430SVaclav Hapla PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[])
6475c6c1daeSBarry Smith {
6485c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
6495c6c1daeSBarry Smith 
6505c6c1daeSBarry Smith   PetscFunctionBegin;
6515c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
652c959eef4SJed Brown   PetscValidPointer(name,2);
653a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
6540298fd71SBarry Smith   else *name = NULL;
6555c6c1daeSBarry Smith   PetscFunctionReturn(0);
6565c6c1daeSBarry Smith }
6575c6c1daeSBarry Smith 
6588c557b5aSVaclav Hapla /*@
659874270d9SVaclav Hapla   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(),
660874270d9SVaclav Hapla   and return this group's ID and file ID.
661874270d9SVaclav Hapla   If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID.
662874270d9SVaclav Hapla 
663874270d9SVaclav Hapla   Not collective
664874270d9SVaclav Hapla 
665874270d9SVaclav Hapla   Input Parameter:
666874270d9SVaclav Hapla . viewer - the PetscViewer
667874270d9SVaclav Hapla 
668874270d9SVaclav Hapla   Output Parameter:
669874270d9SVaclav Hapla + fileId - The HDF5 file ID
670874270d9SVaclav Hapla - groupId - The HDF5 group ID
671874270d9SVaclav Hapla 
672e74616beSVaclav Hapla   Notes:
673e74616beSVaclav Hapla   If the viewer is writable, the group is created if it doesn't exist yet.
674e74616beSVaclav Hapla 
675874270d9SVaclav Hapla   Level: intermediate
676874270d9SVaclav Hapla 
677874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
678874270d9SVaclav Hapla @*/
67954dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
68054dbf706SVaclav Hapla {
68190e03537SVaclav Hapla   hid_t          file_id;
68290e03537SVaclav Hapla   H5O_type_t     type;
683*76d59af2SVaclav Hapla   const char     *groupName = NULL, *fileName = NULL;
684*76d59af2SVaclav Hapla   PetscBool      writable, has;
68554dbf706SVaclav Hapla   PetscErrorCode ierr;
68654dbf706SVaclav Hapla 
68754dbf706SVaclav Hapla   PetscFunctionBegin;
688*76d59af2SVaclav Hapla   ierr = PetscViewerWritable(viewer, &writable);CHKERRQ(ierr);
68954dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr);
690*76d59af2SVaclav Hapla   ierr = PetscViewerFileGetName(viewer, &fileName);CHKERRQ(ierr);
69154dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr);
692*76d59af2SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type);CHKERRQ(ierr);
693*76d59af2SVaclav Hapla   if (!has) {
694*76d59af2SVaclav Hapla     if (!writable) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Group %s does not exist and file %s is not open for writing", groupName, fileName);
695*76d59af2SVaclav Hapla     else           SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName);
696*76d59af2SVaclav Hapla   }
697*76d59af2SVaclav Hapla   if (type != H5O_TYPE_GROUP) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s in file %s resolves to something which is not a group", groupName, fileName);
69890e03537SVaclav Hapla   PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT));
69954dbf706SVaclav Hapla   *fileId  = file_id;
70054dbf706SVaclav Hapla   PetscFunctionReturn(0);
70154dbf706SVaclav Hapla }
70254dbf706SVaclav Hapla 
7033ef9c667SSatish Balay /*@
7045c6c1daeSBarry Smith   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
7055c6c1daeSBarry Smith 
7065c6c1daeSBarry Smith   Not collective
7075c6c1daeSBarry Smith 
7085c6c1daeSBarry Smith   Input Parameter:
7095c6c1daeSBarry Smith . viewer - the PetscViewer
7105c6c1daeSBarry Smith 
7115c6c1daeSBarry Smith   Level: intermediate
7125c6c1daeSBarry Smith 
7135c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
7145c6c1daeSBarry Smith @*/
7155c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
7165c6c1daeSBarry Smith {
7175c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7185c6c1daeSBarry Smith 
7195c6c1daeSBarry Smith   PetscFunctionBegin;
7205c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7215c6c1daeSBarry Smith   ++hdf5->timestep;
7225c6c1daeSBarry Smith   PetscFunctionReturn(0);
7235c6c1daeSBarry Smith }
7245c6c1daeSBarry Smith 
7253ef9c667SSatish Balay /*@
7265c6c1daeSBarry Smith   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
7275c6c1daeSBarry Smith   of -1 disables blocking with timesteps.
7285c6c1daeSBarry Smith 
7295c6c1daeSBarry Smith   Not collective
7305c6c1daeSBarry Smith 
7315c6c1daeSBarry Smith   Input Parameters:
7325c6c1daeSBarry Smith + viewer - the PetscViewer
7335c6c1daeSBarry Smith - timestep - The timestep number
7345c6c1daeSBarry Smith 
7355c6c1daeSBarry Smith   Level: intermediate
7365c6c1daeSBarry Smith 
7375c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
7385c6c1daeSBarry Smith @*/
7395c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
7405c6c1daeSBarry Smith {
7415c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7425c6c1daeSBarry Smith 
7435c6c1daeSBarry Smith   PetscFunctionBegin;
7445c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7455c6c1daeSBarry Smith   hdf5->timestep = timestep;
7465c6c1daeSBarry Smith   PetscFunctionReturn(0);
7475c6c1daeSBarry Smith }
7485c6c1daeSBarry Smith 
7493ef9c667SSatish Balay /*@
7505c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
7515c6c1daeSBarry Smith 
7525c6c1daeSBarry Smith   Not collective
7535c6c1daeSBarry Smith 
7545c6c1daeSBarry Smith   Input Parameter:
7555c6c1daeSBarry Smith . viewer - the PetscViewer
7565c6c1daeSBarry Smith 
7575c6c1daeSBarry Smith   Output Parameter:
7585c6c1daeSBarry Smith . timestep - The timestep number
7595c6c1daeSBarry Smith 
7605c6c1daeSBarry Smith   Level: intermediate
7615c6c1daeSBarry Smith 
7625c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
7635c6c1daeSBarry Smith @*/
7645c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
7655c6c1daeSBarry Smith {
7665c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7675c6c1daeSBarry Smith 
7685c6c1daeSBarry Smith   PetscFunctionBegin;
7695c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7705c6c1daeSBarry Smith   PetscValidPointer(timestep,2);
7715c6c1daeSBarry Smith   *timestep = hdf5->timestep;
7725c6c1daeSBarry Smith   PetscFunctionReturn(0);
7735c6c1daeSBarry Smith }
7745c6c1daeSBarry Smith 
77536bce6e8SMatthew G. Knepley /*@C
77636bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
77736bce6e8SMatthew G. Knepley 
77836bce6e8SMatthew G. Knepley   Not collective
77936bce6e8SMatthew G. Knepley 
78036bce6e8SMatthew G. Knepley   Input Parameter:
78136bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
78236bce6e8SMatthew G. Knepley 
78336bce6e8SMatthew G. Knepley   Output Parameter:
78436bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
78536bce6e8SMatthew G. Knepley 
78636bce6e8SMatthew G. Knepley   Level: advanced
78736bce6e8SMatthew G. Knepley 
78836bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
78936bce6e8SMatthew G. Knepley @*/
79036bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
79136bce6e8SMatthew G. Knepley {
79236bce6e8SMatthew G. Knepley   PetscFunctionBegin;
79336bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
79436bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
79536bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_LLONG;
79636bce6e8SMatthew G. Knepley #else
79736bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_INT;
79836bce6e8SMatthew G. Knepley #endif
79936bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
80036bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
80136bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
802c9450970SBarry Smith   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_INT;
803de3ba262SVaclav Hapla   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_INT;
80436bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
80536bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
80636bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
8077e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
80836bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
80936bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
81036bce6e8SMatthew G. Knepley }
81136bce6e8SMatthew G. Knepley 
81236bce6e8SMatthew G. Knepley /*@C
81336bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
81436bce6e8SMatthew G. Knepley 
81536bce6e8SMatthew G. Knepley   Not collective
81636bce6e8SMatthew G. Knepley 
81736bce6e8SMatthew G. Knepley   Input Parameter:
81836bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
81936bce6e8SMatthew G. Knepley 
82036bce6e8SMatthew G. Knepley   Output Parameter:
82136bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
82236bce6e8SMatthew G. Knepley 
82336bce6e8SMatthew G. Knepley   Level: advanced
82436bce6e8SMatthew G. Knepley 
82536bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
82636bce6e8SMatthew G. Knepley @*/
82736bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
82836bce6e8SMatthew G. Knepley {
82936bce6e8SMatthew G. Knepley   PetscFunctionBegin;
83036bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
83136bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
83236bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
83336bce6e8SMatthew G. Knepley #else
83436bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
83536bce6e8SMatthew G. Knepley #endif
83636bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
83736bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
83836bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
83936bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
84036bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
84136bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
8427e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
84336bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
84436bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
84536bce6e8SMatthew G. Knepley }
84636bce6e8SMatthew G. Knepley 
847df863907SAlex Fikl /*@C
848b8ee85a1SVaclav Hapla  PetscViewerHDF5WriteAttribute - Write an attribute
84936bce6e8SMatthew G. Knepley 
85036bce6e8SMatthew G. Knepley   Input Parameters:
85136bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer
85257d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
85336bce6e8SMatthew G. Knepley . name   - The attribute name
85436bce6e8SMatthew G. Knepley . datatype - The attribute type
85536bce6e8SMatthew G. Knepley - value    - The attribute value
85636bce6e8SMatthew G. Knepley 
85736bce6e8SMatthew G. Knepley   Level: advanced
85836bce6e8SMatthew G. Knepley 
859577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
86036bce6e8SMatthew G. Knepley @*/
86157d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value)
86236bce6e8SMatthew G. Knepley {
86357d22f7dSVaclav Hapla   char           *parent;
86460568592SMatthew G. Knepley   hid_t          h5, dataspace, obj, attribute, dtype;
865080f144cSVaclav Hapla   PetscBool      has;
86636bce6e8SMatthew G. Knepley   PetscErrorCode ierr;
86736bce6e8SMatthew G. Knepley 
86836bce6e8SMatthew G. Knepley   PetscFunctionBegin;
8695cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
87057d22f7dSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset, 2);
871c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
872b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
87357d22f7dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
874bb286ee1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr);
875b67695ecSVaclav Hapla   ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);
87636bce6e8SMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
8777e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
8787e97c476SMatthew G. Knepley     size_t len;
8797e97c476SMatthew G. Knepley     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
880729ab6d9SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
8817e97c476SMatthew G. Knepley   }
88236bce6e8SMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
883729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
88460568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
885080f144cSVaclav Hapla   if (has) {
886080f144cSVaclav Hapla     PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
887080f144cSVaclav Hapla   } else {
88860568592SMatthew G. Knepley     PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
889080f144cSVaclav Hapla   }
890729ab6d9SBarry Smith   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
891729ab6d9SBarry Smith   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
892729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
89360568592SMatthew G. Knepley   PetscStackCallHDF5(H5Oclose,(obj));
894729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(dataspace));
89557d22f7dSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
89636bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
89736bce6e8SMatthew G. Knepley }
89836bce6e8SMatthew G. Knepley 
899df863907SAlex Fikl /*@C
900577e0e04SVaclav Hapla  PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name
901577e0e04SVaclav Hapla 
902577e0e04SVaclav Hapla   Input Parameters:
903577e0e04SVaclav Hapla + viewer   - The HDF5 viewer
904577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
905577e0e04SVaclav Hapla . name     - The attribute name
906577e0e04SVaclav Hapla . datatype - The attribute type
907577e0e04SVaclav Hapla - value    - The attribute value
908577e0e04SVaclav Hapla 
909577e0e04SVaclav Hapla   Notes:
910577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
911577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
912577e0e04SVaclav Hapla 
913577e0e04SVaclav Hapla   Level: advanced
914577e0e04SVaclav Hapla 
915577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
916577e0e04SVaclav Hapla @*/
917577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
918577e0e04SVaclav Hapla {
919577e0e04SVaclav Hapla   PetscErrorCode ierr;
920577e0e04SVaclav Hapla 
921577e0e04SVaclav Hapla   PetscFunctionBegin;
922577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
923577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
924577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
925b7e81f0fSVaclav Hapla   PetscValidPointer(value,5);
926577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
927577e0e04SVaclav Hapla   ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
928577e0e04SVaclav Hapla   PetscFunctionReturn(0);
929577e0e04SVaclav Hapla }
930577e0e04SVaclav Hapla 
931577e0e04SVaclav Hapla /*@C
932b8ee85a1SVaclav Hapla  PetscViewerHDF5ReadAttribute - Read an attribute
933ad0c61feSMatthew G. Knepley 
934ad0c61feSMatthew G. Knepley   Input Parameters:
935ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer
93657d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
937ad0c61feSMatthew G. Knepley . name   - The attribute name
938ad0c61feSMatthew G. Knepley - datatype - The attribute type
939ad0c61feSMatthew G. Knepley 
940ad0c61feSMatthew G. Knepley   Output Parameter:
941ad0c61feSMatthew G. Knepley . value    - The attribute value
942ad0c61feSMatthew G. Knepley 
943708d4cb3SBarry Smith   Notes: If the datatype is PETSC_STRING one must PetscFree() the obtained value when it is no longer needed.
944708d4cb3SBarry Smith 
945ad0c61feSMatthew G. Knepley   Level: advanced
946ad0c61feSMatthew G. Knepley 
947577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
948ad0c61feSMatthew G. Knepley @*/
94957d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value)
950ad0c61feSMatthew G. Knepley {
95157d22f7dSVaclav Hapla   char           *parent;
952f0b58479SMatthew G. Knepley   hid_t          h5, obj, attribute, atype, dtype;
953969748fdSVaclav Hapla   PetscBool      has;
954ad0c61feSMatthew G. Knepley   PetscErrorCode ierr;
955ad0c61feSMatthew G. Knepley 
956ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
9575cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
95857d22f7dSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset, 2);
959c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
960b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
96157d22f7dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
962969748fdSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr);
963969748fdSVaclav Hapla   if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);}
964969748fdSVaclav Hapla   if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name);
965ad0c61feSMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
966ad0c61feSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
96760568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
96860568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
969f0b58479SMatthew G. Knepley   if (datatype == PETSC_STRING) {
970f0b58479SMatthew G. Knepley     size_t len;
971e90c5075SMarek Pecha     PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
97215b861d2SVaclav Hapla     PetscStackCall("H5Tget_size",len = H5Tget_size(atype));
973708d4cb3SBarry Smith     ierr = PetscMalloc((len+1) * sizeof(char), value);CHKERRQ(ierr);
974708d4cb3SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
975708d4cb3SBarry Smith     PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value));
976708d4cb3SBarry Smith   } else {
97770efba33SBarry Smith     PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
978708d4cb3SBarry Smith   }
979729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
980e90c5075SMarek Pecha   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
981e90c5075SMarek Pecha   PetscStackCallHDF5(H5Oclose,(obj));
98257d22f7dSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
983ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
984ad0c61feSMatthew G. Knepley }
985ad0c61feSMatthew G. Knepley 
986577e0e04SVaclav Hapla /*@C
987577e0e04SVaclav Hapla  PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name
988577e0e04SVaclav Hapla 
989577e0e04SVaclav Hapla   Input Parameters:
990577e0e04SVaclav Hapla + viewer   - The HDF5 viewer
991577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
992577e0e04SVaclav Hapla . name     - The attribute name
993577e0e04SVaclav Hapla - datatype - The attribute type
994577e0e04SVaclav Hapla 
995577e0e04SVaclav Hapla   Output Parameter:
996577e0e04SVaclav Hapla . value    - The attribute value
997577e0e04SVaclav Hapla 
998577e0e04SVaclav Hapla   Notes:
999577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1000577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
1001577e0e04SVaclav Hapla 
1002577e0e04SVaclav Hapla   Level: advanced
1003577e0e04SVaclav Hapla 
1004577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1005577e0e04SVaclav Hapla @*/
1006577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value)
1007577e0e04SVaclav Hapla {
1008577e0e04SVaclav Hapla   PetscErrorCode ierr;
1009577e0e04SVaclav Hapla 
1010577e0e04SVaclav Hapla   PetscFunctionBegin;
1011577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1012577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
1013577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
1014b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
1015577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
1016577e0e04SVaclav Hapla   ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
1017577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1018577e0e04SVaclav Hapla }
1019577e0e04SVaclav Hapla 
1020a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
1021a75c3fd4SVaclav Hapla {
1022a75c3fd4SVaclav Hapla   htri_t exists;
1023a75c3fd4SVaclav Hapla   hid_t group;
1024a75c3fd4SVaclav Hapla 
1025a75c3fd4SVaclav Hapla   PetscFunctionBegin;
1026a75c3fd4SVaclav Hapla   PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT));
1027a75c3fd4SVaclav Hapla   if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT));
1028a75c3fd4SVaclav Hapla   if (!exists && createGroup) {
1029a75c3fd4SVaclav Hapla     PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1030a75c3fd4SVaclav Hapla     PetscStackCallHDF5(H5Gclose,(group));
1031a75c3fd4SVaclav Hapla     exists = PETSC_TRUE;
1032a75c3fd4SVaclav Hapla   }
1033a75c3fd4SVaclav Hapla   *exists_ = (PetscBool) exists;
1034a75c3fd4SVaclav Hapla   PetscFunctionReturn(0);
1035a75c3fd4SVaclav Hapla }
1036a75c3fd4SVaclav Hapla 
1037bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
10385cdeb108SMatthew G. Knepley {
103990e03537SVaclav Hapla   const char     rootGroupName[] = "/";
10405cdeb108SMatthew G. Knepley   hid_t          h5;
1041e5bf9ebcSVaclav Hapla   PetscBool      exists=PETSC_FALSE;
104215b861d2SVaclav Hapla   PetscInt       i;
104315b861d2SVaclav Hapla   int            n;
104485688ae2SVaclav Hapla   char           **hierarchy;
104585688ae2SVaclav Hapla   char           buf[PETSC_MAX_PATH_LEN]="";
10465cdeb108SMatthew G. Knepley   PetscErrorCode ierr;
10475cdeb108SMatthew G. Knepley 
10485cdeb108SMatthew G. Knepley   PetscFunctionBegin;
10495cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
105090e03537SVaclav Hapla   if (name) PetscValidCharPointer(name, 2);
105190e03537SVaclav Hapla   else name = rootGroupName;
1052ccd66a83SVaclav Hapla   if (has) {
105356cc0592SVaclav Hapla     PetscValidIntPointer(has, 3);
10545cdeb108SMatthew G. Knepley     *has = PETSC_FALSE;
1055ccd66a83SVaclav Hapla   }
1056ccd66a83SVaclav Hapla   if (otype) {
1057ccd66a83SVaclav Hapla     PetscValidIntPointer(otype, 4);
105856cc0592SVaclav Hapla     *otype = H5O_TYPE_UNKNOWN;
1059ccd66a83SVaclav Hapla   }
1060c0cd0301SJed Brown   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
106185688ae2SVaclav Hapla 
106285688ae2SVaclav Hapla   /*
106385688ae2SVaclav Hapla      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
106485688ae2SVaclav Hapla      Hence, each of them needs to be tested separately:
106585688ae2SVaclav Hapla      1) whether it's a valid link
106685688ae2SVaclav Hapla      2) whether this link resolves to an object
106785688ae2SVaclav Hapla      See H5Oexists_by_name() documentation.
106885688ae2SVaclav Hapla   */
106985688ae2SVaclav Hapla   ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr);
107085688ae2SVaclav Hapla   if (!n) {
107185688ae2SVaclav Hapla     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1072ccd66a83SVaclav Hapla     if (has)   *has   = PETSC_TRUE;
1073ccd66a83SVaclav Hapla     if (otype) *otype = H5O_TYPE_GROUP;
107485688ae2SVaclav Hapla     ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
107585688ae2SVaclav Hapla     PetscFunctionReturn(0);
107685688ae2SVaclav Hapla   }
107785688ae2SVaclav Hapla   for (i=0; i<n; i++) {
107885688ae2SVaclav Hapla     ierr = PetscStrcat(buf,"/");CHKERRQ(ierr);
107985688ae2SVaclav Hapla     ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr);
1080a75c3fd4SVaclav Hapla     ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr);
1081a75c3fd4SVaclav Hapla     if (!exists) break;
108285688ae2SVaclav Hapla   }
108385688ae2SVaclav Hapla   ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
108485688ae2SVaclav Hapla 
108585688ae2SVaclav Hapla   /* If the object exists, get its type */
1086ccd66a83SVaclav Hapla   if (exists && otype) {
10875cdeb108SMatthew G. Knepley     H5O_info_t info;
10885cdeb108SMatthew G. Knepley 
108976276748SVaclav Hapla     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
109004633f7fSVaclav Hapla     PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT));
109156cc0592SVaclav Hapla     *otype = info.type;
10925cdeb108SMatthew G. Knepley   }
1093ccd66a83SVaclav Hapla   if (has) *has = exists;
10945cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
10955cdeb108SMatthew G. Knepley }
10965cdeb108SMatthew G. Knepley 
1097ecce1506SVaclav Hapla /*@
10980a9f38efSVaclav Hapla  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
10990a9f38efSVaclav Hapla 
11000a9f38efSVaclav Hapla   Input Parameters:
11010a9f38efSVaclav Hapla . viewer - The HDF5 viewer
11020a9f38efSVaclav Hapla 
11030a9f38efSVaclav Hapla   Output Parameter:
11040a9f38efSVaclav Hapla . has    - Flag for group existence
11050a9f38efSVaclav Hapla 
11060a9f38efSVaclav Hapla   Notes:
11070a9f38efSVaclav Hapla   If the path exists but is not a group, this returns PETSC_FALSE as well.
11080a9f38efSVaclav Hapla 
11090a9f38efSVaclav Hapla   Level: advanced
11100a9f38efSVaclav Hapla 
11110a9f38efSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup()
11120a9f38efSVaclav Hapla @*/
11130a9f38efSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has)
11140a9f38efSVaclav Hapla {
11150a9f38efSVaclav Hapla   H5O_type_t type;
11160a9f38efSVaclav Hapla   const char *name;
11170a9f38efSVaclav Hapla   PetscErrorCode ierr;
11180a9f38efSVaclav Hapla 
11190a9f38efSVaclav Hapla   PetscFunctionBegin;
11200a9f38efSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
11210a9f38efSVaclav Hapla   PetscValidIntPointer(has,2);
11220a9f38efSVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr);
11230a9f38efSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr);
11240a9f38efSVaclav Hapla   *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE;
11250a9f38efSVaclav Hapla   PetscFunctionReturn(0);
11260a9f38efSVaclav Hapla }
11270a9f38efSVaclav Hapla 
11280a9f38efSVaclav Hapla /*@
1129e3f143f7SVaclav Hapla  PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1130ecce1506SVaclav Hapla 
1131ecce1506SVaclav Hapla   Input Parameters:
1132ecce1506SVaclav Hapla + viewer - The HDF5 viewer
1133ecce1506SVaclav Hapla - obj    - The named object
1134ecce1506SVaclav Hapla 
1135ecce1506SVaclav Hapla   Output Parameter:
1136ecce1506SVaclav Hapla . has    - Flag for dataset existence; PETSC_FALSE for unnamed object
1137ecce1506SVaclav Hapla 
1138e3f143f7SVaclav Hapla   Notes:
1139e3f143f7SVaclav Hapla   If the path exists but is not a dataset, this returns PETSC_FALSE as well.
1140e3f143f7SVaclav Hapla 
1141ecce1506SVaclav Hapla   Level: advanced
1142ecce1506SVaclav Hapla 
1143e3f143f7SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1144ecce1506SVaclav Hapla @*/
1145ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1146ecce1506SVaclav Hapla {
114756cc0592SVaclav Hapla   H5O_type_t type;
11486c132bc1SVaclav Hapla   char *path;
1149ecce1506SVaclav Hapla   PetscErrorCode ierr;
1150ecce1506SVaclav Hapla 
1151ecce1506SVaclav Hapla   PetscFunctionBegin;
1152c57a1a9aSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1153c57a1a9aSVaclav Hapla   PetscValidHeader(obj,2);
1154c57a1a9aSVaclav Hapla   PetscValidIntPointer(has,3);
1155ecce1506SVaclav Hapla   *has = PETSC_FALSE;
1156e3f143f7SVaclav Hapla   if (!obj->name) PetscFunctionReturn(0);
11576c132bc1SVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr);
11586c132bc1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr);
115956cc0592SVaclav Hapla   *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE;
11606c132bc1SVaclav Hapla   ierr = PetscFree(path);CHKERRQ(ierr);
1161ecce1506SVaclav Hapla   PetscFunctionReturn(0);
1162ecce1506SVaclav Hapla }
1163ecce1506SVaclav Hapla 
1164df863907SAlex Fikl /*@C
1165ece83b6cSVaclav Hapla  PetscViewerHDF5HasAttribute - Check whether an attribute exists
1166e7bdbf8eSMatthew G. Knepley 
1167e7bdbf8eSMatthew G. Knepley   Input Parameters:
1168e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer
116910e69e8fSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
1170e7bdbf8eSMatthew G. Knepley - name   - The attribute name
1171e7bdbf8eSMatthew G. Knepley 
1172e7bdbf8eSMatthew G. Knepley   Output Parameter:
1173e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
1174e7bdbf8eSMatthew G. Knepley 
1175e7bdbf8eSMatthew G. Knepley   Level: advanced
1176e7bdbf8eSMatthew G. Knepley 
1177577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1178e7bdbf8eSMatthew G. Knepley @*/
117910e69e8fSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has)
1180e7bdbf8eSMatthew G. Knepley {
118110e69e8fSVaclav Hapla   char           *parent;
1182e7bdbf8eSMatthew G. Knepley   PetscErrorCode ierr;
1183e7bdbf8eSMatthew G. Knepley 
1184e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
11855cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
118610e69e8fSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset,2);
1187c57a1a9aSVaclav Hapla   PetscValidCharPointer(name,3);
1188c57a1a9aSVaclav Hapla   PetscValidIntPointer(has,4);
118910e69e8fSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
1190bb286ee1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr);
119110e69e8fSVaclav Hapla   if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);}
119210e69e8fSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
119306db490cSVaclav Hapla   PetscFunctionReturn(0);
119406db490cSVaclav Hapla }
119506db490cSVaclav Hapla 
1196577e0e04SVaclav Hapla /*@C
1197577e0e04SVaclav Hapla  PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name
1198577e0e04SVaclav Hapla 
1199577e0e04SVaclav Hapla   Input Parameters:
1200577e0e04SVaclav Hapla + viewer - The HDF5 viewer
1201577e0e04SVaclav Hapla . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1202577e0e04SVaclav Hapla - name   - The attribute name
1203577e0e04SVaclav Hapla 
1204577e0e04SVaclav Hapla   Output Parameter:
1205577e0e04SVaclav Hapla . has    - Flag for attribute existence
1206577e0e04SVaclav Hapla 
1207577e0e04SVaclav Hapla   Notes:
1208577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1209577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
1210577e0e04SVaclav Hapla 
1211577e0e04SVaclav Hapla   Level: advanced
1212577e0e04SVaclav Hapla 
1213577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1214577e0e04SVaclav Hapla @*/
1215577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1216577e0e04SVaclav Hapla {
1217577e0e04SVaclav Hapla   PetscErrorCode ierr;
1218577e0e04SVaclav Hapla 
1219577e0e04SVaclav Hapla   PetscFunctionBegin;
1220577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1221577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
1222577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
1223577e0e04SVaclav Hapla   PetscValidIntPointer(has,4);
1224577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
1225577e0e04SVaclav Hapla   ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr);
1226577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1227577e0e04SVaclav Hapla }
1228577e0e04SVaclav Hapla 
122906db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
123006db490cSVaclav Hapla {
123106db490cSVaclav Hapla   hid_t          h5;
123206db490cSVaclav Hapla   htri_t         hhas;
123306db490cSVaclav Hapla   PetscErrorCode ierr;
123406db490cSVaclav Hapla 
123506db490cSVaclav Hapla   PetscFunctionBegin;
123606db490cSVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
12372f936e54SVaclav Hapla   PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT));
12385cdeb108SMatthew G. Knepley   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1239e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
1240e7bdbf8eSMatthew G. Knepley }
1241e7bdbf8eSMatthew G. Knepley 
1242a75e6a4aSMatthew G. Knepley /*
1243a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1244a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
1245a75e6a4aSMatthew G. Knepley */
1246d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1247a75e6a4aSMatthew G. Knepley 
1248a75e6a4aSMatthew G. Knepley /*@C
1249a75e6a4aSMatthew G. Knepley   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
1250a75e6a4aSMatthew G. Knepley 
1251d083f849SBarry Smith   Collective
1252a75e6a4aSMatthew G. Knepley 
1253a75e6a4aSMatthew G. Knepley   Input Parameter:
1254a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer
1255a75e6a4aSMatthew G. Knepley 
1256a75e6a4aSMatthew G. Knepley   Level: intermediate
1257a75e6a4aSMatthew G. Knepley 
1258a75e6a4aSMatthew G. Knepley   Options Database Keys:
125910699b91SBarry Smith . -viewer_hdf5_filename <name> - name of the HDF5 file
1260a75e6a4aSMatthew G. Knepley 
1261a75e6a4aSMatthew G. Knepley   Environmental variables:
126210699b91SBarry Smith . PETSC_VIEWER_HDF5_FILENAME - name of the HDF5 file
1263a75e6a4aSMatthew G. Knepley 
1264a75e6a4aSMatthew G. Knepley   Notes:
1265a75e6a4aSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
1266a75e6a4aSMatthew G. Knepley   an error code.  The HDF5 PetscViewer is usually used in the form
1267a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1268a75e6a4aSMatthew G. Knepley 
1269a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1270a75e6a4aSMatthew G. Knepley @*/
1271a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1272a75e6a4aSMatthew G. Knepley {
1273a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
1274a75e6a4aSMatthew G. Knepley   PetscBool      flg;
1275a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
1276a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
1277a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
1278a75e6a4aSMatthew G. Knepley 
1279a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
1280908793a3SLisandro 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);}
1281a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1282908793a3SLisandro Dalcin     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,NULL);
1283908793a3SLisandro Dalcin     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1284a75e6a4aSMatthew G. Knepley   }
128547435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1286908793a3SLisandro Dalcin   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1287a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
1288a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
12892cb5e1ccSBarry Smith     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1290a75e6a4aSMatthew G. Knepley     if (!flg) {
1291a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname,"output.h5");
12922cb5e1ccSBarry Smith       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1293a75e6a4aSMatthew G. Knepley     }
1294a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
12952cb5e1ccSBarry Smith     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1296a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
12972cb5e1ccSBarry Smith     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
129847435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
1299908793a3SLisandro Dalcin     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1300a75e6a4aSMatthew G. Knepley   }
1301a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
13022cb5e1ccSBarry Smith   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1303a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
1304a75e6a4aSMatthew G. Knepley }
1305