xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 0b62783d0d5822e7db3bde533443071fe8ef06fc)
1af0996ceSBarry Smith #include <petsc/private/viewerimpl.h>    /*I   "petscsys.h"   I*/
2d70abbfaSBarry Smith #include <petscviewerhdf5.h>    /*I   "petscviewerhdf5.h"   I*/
3ded06c3fSVaclav Hapla #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE < 10800)
4ded06c3fSVaclav Hapla #error "PETSc needs HDF5 version >= 1.8.0"
5ded06c3fSVaclav Hapla #endif
65c6c1daeSBarry Smith 
7bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool*, H5O_type_t*);
806db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool*);
906db490cSVaclav Hapla 
105c6c1daeSBarry Smith typedef struct GroupList {
115c6c1daeSBarry Smith   const char       *name;
125c6c1daeSBarry Smith   struct GroupList *next;
135c6c1daeSBarry Smith } GroupList;
145c6c1daeSBarry Smith 
155c6c1daeSBarry Smith typedef struct {
165c6c1daeSBarry Smith   char          *filename;
175c6c1daeSBarry Smith   PetscFileMode btype;
185c6c1daeSBarry Smith   hid_t         file_id;
195c6c1daeSBarry Smith   PetscInt      timestep;
205c6c1daeSBarry Smith   GroupList     *groups;
2182ea9b62SBarry Smith   PetscBool     basedimension2;  /* save vectors and DMDA vectors with a dimension of at least 2 even if the bs/dof is 1 */
229a0502c6SHåkon Strandenes   PetscBool     spoutput;  /* write data in single precision even if PETSc is compiled with double precision PetscReal */
23058bd781SVaclav Hapla   char          *mataij_iname;
24058bd781SVaclav Hapla   char          *mataij_jname;
25058bd781SVaclav Hapla   char          *mataij_aname;
26058bd781SVaclav Hapla   char          *mataij_cname;
27b723ab35SVaclav Hapla   PetscBool     mataij_names_set;
285c6c1daeSBarry Smith } PetscViewer_HDF5;
295c6c1daeSBarry Smith 
30eb5a92b4SVaclav Hapla struct _n_HDF5ReadCtx {
31eb5a92b4SVaclav Hapla   hid_t file, group, dataset, dataspace, plist;
32eb5a92b4SVaclav Hapla   PetscInt timestep;
3309dabeb0SVaclav Hapla   PetscBool complexVal, dim2, horizontal;
34eb5a92b4SVaclav Hapla };
35eb5a92b4SVaclav Hapla typedef struct _n_HDF5ReadCtx* HDF5ReadCtx;
36eb5a92b4SVaclav Hapla 
376c132bc1SVaclav Hapla static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char objname[], char **fullpath)
386c132bc1SVaclav Hapla {
396c132bc1SVaclav Hapla   const char *group;
406c132bc1SVaclav Hapla   char buf[PETSC_MAX_PATH_LEN]="";
416c132bc1SVaclav Hapla   PetscErrorCode ierr;
426c132bc1SVaclav Hapla 
436c132bc1SVaclav Hapla   PetscFunctionBegin;
446c132bc1SVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr);
456c132bc1SVaclav Hapla   ierr = PetscStrcat(buf, group);CHKERRQ(ierr);
466c132bc1SVaclav Hapla   ierr = PetscStrcat(buf, "/");CHKERRQ(ierr);
476c132bc1SVaclav Hapla   ierr = PetscStrcat(buf, objname);CHKERRQ(ierr);
486c132bc1SVaclav Hapla   ierr = PetscStrallocpy(buf, fullpath);CHKERRQ(ierr);
496c132bc1SVaclav Hapla   PetscFunctionReturn(0);
506c132bc1SVaclav Hapla }
516c132bc1SVaclav Hapla 
52577e0e04SVaclav Hapla static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj)
53577e0e04SVaclav Hapla {
54577e0e04SVaclav Hapla   PetscBool has;
55577e0e04SVaclav Hapla   const char *group;
56577e0e04SVaclav Hapla   PetscErrorCode ierr;
57577e0e04SVaclav Hapla 
58577e0e04SVaclav Hapla   PetscFunctionBegin;
59577e0e04SVaclav Hapla   if (!obj->name) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
60577e0e04SVaclav Hapla   ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr);
61577e0e04SVaclav Hapla   if (!has) {
62577e0e04SVaclav Hapla     ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr);
63577e0e04SVaclav Hapla     SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) %s not stored in group %s", obj->name, group);
64577e0e04SVaclav Hapla   }
65577e0e04SVaclav Hapla   PetscFunctionReturn(0);
66577e0e04SVaclav Hapla }
67577e0e04SVaclav Hapla 
684416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v)
6982ea9b62SBarry Smith {
7082ea9b62SBarry Smith   PetscErrorCode   ierr;
7182ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;
7282ea9b62SBarry Smith 
7382ea9b62SBarry Smith   PetscFunctionBegin;
7482ea9b62SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr);
7582ea9b62SBarry Smith   ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr);
769a0502c6SHåkon Strandenes   ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr);
7782ea9b62SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
7882ea9b62SBarry Smith   PetscFunctionReturn(0);
7982ea9b62SBarry Smith }
8082ea9b62SBarry Smith 
815c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
825c6c1daeSBarry Smith {
835c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
845c6c1daeSBarry Smith   PetscErrorCode   ierr;
855c6c1daeSBarry Smith 
865c6c1daeSBarry Smith   PetscFunctionBegin;
875c6c1daeSBarry Smith   ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);
88729ab6d9SBarry Smith   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
895c6c1daeSBarry Smith   PetscFunctionReturn(0);
905c6c1daeSBarry Smith }
915c6c1daeSBarry Smith 
929b1403beSVaclav Hapla static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
935c6c1daeSBarry Smith {
945c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
955c6c1daeSBarry Smith   PetscErrorCode   ierr;
965c6c1daeSBarry Smith 
975c6c1daeSBarry Smith   PetscFunctionBegin;
985c6c1daeSBarry Smith   ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr);
995c6c1daeSBarry Smith   while (hdf5->groups) {
1005c6c1daeSBarry Smith     GroupList *tmp = hdf5->groups->next;
1015c6c1daeSBarry Smith 
1025c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups->name);CHKERRQ(ierr);
1035c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups);CHKERRQ(ierr);
1045c6c1daeSBarry Smith     hdf5->groups = tmp;
1055c6c1daeSBarry Smith   }
10661912a0cSVaclav Hapla   ierr = PetscFree(hdf5->mataij_iname);CHKERRQ(ierr);
10761912a0cSVaclav Hapla   ierr = PetscFree(hdf5->mataij_jname);CHKERRQ(ierr);
10861912a0cSVaclav Hapla   ierr = PetscFree(hdf5->mataij_aname);CHKERRQ(ierr);
10961912a0cSVaclav Hapla   ierr = PetscFree(hdf5->mataij_cname);CHKERRQ(ierr);
1105c6c1daeSBarry Smith   ierr = PetscFree(hdf5);CHKERRQ(ierr);
1110b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
112d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr);
1130b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
114058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr);
115058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr);
116058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetAIJNames_C",NULL);CHKERRQ(ierr);
117058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5GetAIJNames_C",NULL);CHKERRQ(ierr);
1185c6c1daeSBarry Smith   PetscFunctionReturn(0);
1195c6c1daeSBarry Smith }
1205c6c1daeSBarry Smith 
1219b1403beSVaclav Hapla static PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
1225c6c1daeSBarry Smith {
1235c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1245c6c1daeSBarry Smith 
1255c6c1daeSBarry Smith   PetscFunctionBegin;
1265c6c1daeSBarry Smith   hdf5->btype = type;
1275c6c1daeSBarry Smith   PetscFunctionReturn(0);
1285c6c1daeSBarry Smith }
1295c6c1daeSBarry Smith 
130*0b62783dSVaclav Hapla static PetscErrorCode  PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
131*0b62783dSVaclav Hapla {
132*0b62783dSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
133*0b62783dSVaclav Hapla 
134*0b62783dSVaclav Hapla   PetscFunctionBegin;
135*0b62783dSVaclav Hapla   *type = hdf5->btype;
136*0b62783dSVaclav Hapla   PetscFunctionReturn(0);
137*0b62783dSVaclav Hapla }
138*0b62783dSVaclav Hapla 
1399b1403beSVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
14082ea9b62SBarry Smith {
14182ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
14282ea9b62SBarry Smith 
14382ea9b62SBarry Smith   PetscFunctionBegin;
14482ea9b62SBarry Smith   hdf5->basedimension2 = flg;
14582ea9b62SBarry Smith   PetscFunctionReturn(0);
14682ea9b62SBarry Smith }
14782ea9b62SBarry Smith 
148df863907SAlex Fikl /*@
14982ea9b62SBarry Smith      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
15082ea9b62SBarry Smith        dimension of 2.
15182ea9b62SBarry Smith 
15282ea9b62SBarry Smith     Logically Collective on PetscViewer
15382ea9b62SBarry Smith 
15482ea9b62SBarry Smith   Input Parameters:
15582ea9b62SBarry Smith +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
15682ea9b62SBarry 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
15782ea9b62SBarry Smith 
15882ea9b62SBarry Smith   Options Database:
15982ea9b62SBarry 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
16082ea9b62SBarry Smith 
16182ea9b62SBarry Smith 
16295452b02SPatrick Sanan   Notes:
16395452b02SPatrick 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
16482ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
16582ea9b62SBarry Smith 
16682ea9b62SBarry Smith   Level: intermediate
16782ea9b62SBarry Smith 
16882ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
16982ea9b62SBarry Smith 
17082ea9b62SBarry Smith @*/
17182ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
17282ea9b62SBarry Smith {
17382ea9b62SBarry Smith   PetscErrorCode ierr;
17482ea9b62SBarry Smith 
17582ea9b62SBarry Smith   PetscFunctionBegin;
17682ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
17782ea9b62SBarry Smith   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
17882ea9b62SBarry Smith   PetscFunctionReturn(0);
17982ea9b62SBarry Smith }
18082ea9b62SBarry Smith 
181df863907SAlex Fikl /*@
18282ea9b62SBarry Smith      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
18382ea9b62SBarry Smith        dimension of 2.
18482ea9b62SBarry Smith 
18582ea9b62SBarry Smith     Logically Collective on PetscViewer
18682ea9b62SBarry Smith 
18782ea9b62SBarry Smith   Input Parameter:
18882ea9b62SBarry Smith .  viewer - the PetscViewer, must be of type HDF5
18982ea9b62SBarry Smith 
19082ea9b62SBarry Smith   Output Parameter:
19182ea9b62SBarry 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
19282ea9b62SBarry Smith 
19395452b02SPatrick Sanan   Notes:
19495452b02SPatrick 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
19582ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
19682ea9b62SBarry Smith 
19782ea9b62SBarry Smith   Level: intermediate
19882ea9b62SBarry Smith 
19982ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
20082ea9b62SBarry Smith 
20182ea9b62SBarry Smith @*/
20282ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
20382ea9b62SBarry Smith {
20482ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
20582ea9b62SBarry Smith 
20682ea9b62SBarry Smith   PetscFunctionBegin;
20782ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
20882ea9b62SBarry Smith   *flg = hdf5->basedimension2;
20982ea9b62SBarry Smith   PetscFunctionReturn(0);
21082ea9b62SBarry Smith }
21182ea9b62SBarry Smith 
2129b1403beSVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
2139a0502c6SHåkon Strandenes {
2149a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2159a0502c6SHåkon Strandenes 
2169a0502c6SHåkon Strandenes   PetscFunctionBegin;
2179a0502c6SHåkon Strandenes   hdf5->spoutput = flg;
2189a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2199a0502c6SHåkon Strandenes }
2209a0502c6SHåkon Strandenes 
221df863907SAlex Fikl /*@
2229a0502c6SHåkon Strandenes      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
2239a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2249a0502c6SHåkon Strandenes 
2259a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2269a0502c6SHåkon Strandenes 
2279a0502c6SHåkon Strandenes   Input Parameters:
2289a0502c6SHåkon Strandenes +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
2299a0502c6SHåkon Strandenes -  flg - if PETSC_TRUE the data will be written to disk with single precision
2309a0502c6SHåkon Strandenes 
2319a0502c6SHåkon Strandenes   Options Database:
2329a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
2339a0502c6SHåkon Strandenes 
2349a0502c6SHåkon Strandenes 
23595452b02SPatrick Sanan   Notes:
23695452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2379a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2389a0502c6SHåkon Strandenes 
2399a0502c6SHåkon Strandenes   Level: intermediate
2409a0502c6SHåkon Strandenes 
2419a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2429a0502c6SHåkon Strandenes           PetscReal
2439a0502c6SHåkon Strandenes 
2449a0502c6SHåkon Strandenes @*/
2459a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
2469a0502c6SHåkon Strandenes {
2479a0502c6SHåkon Strandenes   PetscErrorCode ierr;
2489a0502c6SHåkon Strandenes 
2499a0502c6SHåkon Strandenes   PetscFunctionBegin;
2509a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2519a0502c6SHåkon Strandenes   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
2529a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2539a0502c6SHåkon Strandenes }
2549a0502c6SHåkon Strandenes 
255df863907SAlex Fikl /*@
2569a0502c6SHåkon Strandenes      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
2579a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2589a0502c6SHåkon Strandenes 
2599a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2609a0502c6SHåkon Strandenes 
2619a0502c6SHåkon Strandenes   Input Parameter:
2629a0502c6SHåkon Strandenes .  viewer - the PetscViewer, must be of type HDF5
2639a0502c6SHåkon Strandenes 
2649a0502c6SHåkon Strandenes   Output Parameter:
2659a0502c6SHåkon Strandenes .  flg - if PETSC_TRUE the data will be written to disk with single precision
2669a0502c6SHåkon Strandenes 
26795452b02SPatrick Sanan   Notes:
26895452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2699a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2709a0502c6SHåkon Strandenes 
2719a0502c6SHåkon Strandenes   Level: intermediate
2729a0502c6SHåkon Strandenes 
2739a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2749a0502c6SHåkon Strandenes           PetscReal
2759a0502c6SHåkon Strandenes 
2769a0502c6SHåkon Strandenes @*/
2779a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
2789a0502c6SHåkon Strandenes {
2799a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2809a0502c6SHåkon Strandenes 
2819a0502c6SHåkon Strandenes   PetscFunctionBegin;
2829a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2839a0502c6SHåkon Strandenes   *flg = hdf5->spoutput;
2849a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2859a0502c6SHåkon Strandenes }
2869a0502c6SHåkon Strandenes 
2879b1403beSVaclav Hapla static PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
2885c6c1daeSBarry Smith {
2895c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2905c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
2915c6c1daeSBarry Smith   MPI_Info          info = MPI_INFO_NULL;
2925c6c1daeSBarry Smith #endif
2935c6c1daeSBarry Smith   hid_t             plist_id;
2945c6c1daeSBarry Smith   PetscErrorCode    ierr;
2955c6c1daeSBarry Smith 
2965c6c1daeSBarry Smith   PetscFunctionBegin;
297f683cc58SToby Isaac   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
298f683cc58SToby Isaac   if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);}
2995c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
3005c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
301729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
3025c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
303729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info));
3045c6c1daeSBarry Smith #endif
3055c6c1daeSBarry Smith   /* Create or open the file collectively */
3065c6c1daeSBarry Smith   switch (hdf5->btype) {
3075c6c1daeSBarry Smith   case FILE_MODE_READ:
308729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
3095c6c1daeSBarry Smith     break;
3105c6c1daeSBarry Smith   case FILE_MODE_APPEND:
311729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
3125c6c1daeSBarry Smith     break;
3135c6c1daeSBarry Smith   case FILE_MODE_WRITE:
314729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
3155c6c1daeSBarry Smith     break;
3165c6c1daeSBarry Smith   default:
3175c6c1daeSBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
3185c6c1daeSBarry Smith   }
3195c6c1daeSBarry Smith   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
320729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
3215c6c1daeSBarry Smith   PetscFunctionReturn(0);
3225c6c1daeSBarry Smith }
3235c6c1daeSBarry Smith 
324d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
325d1232d7fSToby Isaac {
326d1232d7fSToby Isaac   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
327d1232d7fSToby Isaac 
328d1232d7fSToby Isaac   PetscFunctionBegin;
329d1232d7fSToby Isaac   *name = vhdf5->filename;
330d1232d7fSToby Isaac   PetscFunctionReturn(0);
331d1232d7fSToby Isaac }
332d1232d7fSToby Isaac 
3339b1403beSVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetAIJNames_HDF5(PetscViewer viewer, const char iname[], const char jname[], const char aname[], const char cname[])
334058bd781SVaclav Hapla {
335058bd781SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
336058bd781SVaclav Hapla   PetscErrorCode ierr;
337058bd781SVaclav Hapla 
338058bd781SVaclav Hapla   PetscFunctionBegin;
339058bd781SVaclav Hapla   ierr = PetscFree(hdf5->mataij_iname);CHKERRQ(ierr);
340058bd781SVaclav Hapla   ierr = PetscFree(hdf5->mataij_jname);CHKERRQ(ierr);
341058bd781SVaclav Hapla   ierr = PetscFree(hdf5->mataij_aname);CHKERRQ(ierr);
342058bd781SVaclav Hapla   ierr = PetscFree(hdf5->mataij_cname);CHKERRQ(ierr);
343058bd781SVaclav Hapla   ierr = PetscStrallocpy(iname,&hdf5->mataij_iname);CHKERRQ(ierr);
344058bd781SVaclav Hapla   ierr = PetscStrallocpy(jname,&hdf5->mataij_jname);CHKERRQ(ierr);
345058bd781SVaclav Hapla   ierr = PetscStrallocpy(aname,&hdf5->mataij_aname);CHKERRQ(ierr);
346058bd781SVaclav Hapla   ierr = PetscStrallocpy(cname,&hdf5->mataij_cname);CHKERRQ(ierr);
347b723ab35SVaclav Hapla   hdf5->mataij_names_set = PETSC_TRUE;
348058bd781SVaclav Hapla   PetscFunctionReturn(0);
349058bd781SVaclav Hapla }
350058bd781SVaclav Hapla 
351058bd781SVaclav Hapla /*@C
352058bd781SVaclav Hapla   PetscViewerHDF5SetAIJNames - Set the names of the datasets representing the three AIJ (CRS) arrays and the name of the attribute storing the number of columns within the HDF5 file.
353058bd781SVaclav Hapla 
354058bd781SVaclav Hapla   Collective on PetscViewer
355058bd781SVaclav Hapla 
356058bd781SVaclav Hapla   Input Parameters:
357058bd781SVaclav Hapla +  viewer - the PetscViewer; either ASCII or binary
358058bd781SVaclav Hapla .  iname - name of dataset i representing row pointers; that is i[0] = 0, i[row] = i[row-1] + number of elements in that row of the matrix
359058bd781SVaclav Hapla .  jname - name of dataset j representing column indices
360058bd781SVaclav Hapla .  aname - name of dataset a representing matrix values
361058bd781SVaclav Hapla -  cname - name of attribute stoting column count
362058bd781SVaclav Hapla 
363058bd781SVaclav Hapla   Level: advanced
364058bd781SVaclav Hapla 
365058bd781SVaclav Hapla   Notes:
366cbb4c999SVaclav Hapla   Current defaults are (iname, jname, aname, cname) = ("i", "j", "a", "ncols").
367cbb4c999SVaclav Hapla   For PetscViewerFormat PETSC_VIEWER_HDF5_MAT they are ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be loaded.
368058bd781SVaclav Hapla 
369058bd781SVaclav Hapla .seealso: MatLoad(), PetscViewerCreate(), PetscViewerSetType(), PETSCVIEWERHDF5, PetscViewerHDF5GetAIJNames()
370058bd781SVaclav Hapla @*/
371058bd781SVaclav Hapla PetscErrorCode  PetscViewerHDF5SetAIJNames(PetscViewer viewer, const char iname[], const char jname[], const char aname[], const char cname[])
372058bd781SVaclav Hapla {
373058bd781SVaclav Hapla   PetscErrorCode ierr;
374058bd781SVaclav Hapla 
375058bd781SVaclav Hapla   PetscFunctionBegin;
376058bd781SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
377058bd781SVaclav Hapla   PetscValidCharPointer(iname,2);
378058bd781SVaclav Hapla   PetscValidCharPointer(jname,3);
379058bd781SVaclav Hapla   PetscValidCharPointer(aname,4);
380058bd781SVaclav Hapla   PetscValidCharPointer(cname,5);
381058bd781SVaclav Hapla   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetAIJNames_C",(PetscViewer,const char[],const char[],const char[],const char[]),(viewer,iname,jname,aname,cname));CHKERRQ(ierr);
382058bd781SVaclav Hapla   PetscFunctionReturn(0);
383058bd781SVaclav Hapla }
384058bd781SVaclav Hapla 
3859b1403beSVaclav Hapla static PetscErrorCode  PetscViewerHDF5GetAIJNames_HDF5(PetscViewer viewer, const char *iname[], const char *jname[], const char *aname[], const char *cname[])
386058bd781SVaclav Hapla {
387058bd781SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
388058bd781SVaclav Hapla 
389058bd781SVaclav Hapla   PetscFunctionBegin;
390058bd781SVaclav Hapla   *iname = hdf5->mataij_iname;
391058bd781SVaclav Hapla   *jname = hdf5->mataij_jname;
392058bd781SVaclav Hapla   *aname = hdf5->mataij_aname;
393058bd781SVaclav Hapla   *cname = hdf5->mataij_cname;
394058bd781SVaclav Hapla   PetscFunctionReturn(0);
395058bd781SVaclav Hapla }
396058bd781SVaclav Hapla 
397058bd781SVaclav Hapla /*@C
398058bd781SVaclav Hapla   PetscViewerHDF5GetAIJNames - Get the names of the datasets representing the three AIJ (CRS) arrays and the name of the attribute storing the number of columns within the HDF5 file.
399058bd781SVaclav Hapla 
400058bd781SVaclav Hapla   Collective on PetscViewer
401058bd781SVaclav Hapla 
402058bd781SVaclav Hapla   Input Parameters:
403058bd781SVaclav Hapla .  viewer - the PetscViewer; either ASCII or binary
404058bd781SVaclav Hapla 
405058bd781SVaclav Hapla   Output Parameters:
406058bd781SVaclav Hapla +  iname - name of dataset i representing row pointers; that is i[0] = 0, i[row] = i[row-1] + number of elements in that row of the matrix
407058bd781SVaclav Hapla .  jname - name of dataset j representing column indices
408058bd781SVaclav Hapla .  aname - name of dataset a representing matrix values
409058bd781SVaclav Hapla -  cname - name of attribute stoting column count
410058bd781SVaclav Hapla 
411058bd781SVaclav Hapla   Level: advanced
412058bd781SVaclav Hapla 
413058bd781SVaclav Hapla   Notes:
414cbb4c999SVaclav Hapla   Current defaults are (iname, jname, aname, cname) = ("i", "j", "a", "ncols").
415cbb4c999SVaclav Hapla   For PetscViewerFormat PETSC_VIEWER_HDF5_MAT they are ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be loaded.
416058bd781SVaclav Hapla 
417058bd781SVaclav Hapla .seealso: MatLoad(), PetscViewerCreate(), PetscViewerSetType(), PETSCVIEWERHDF5, PetscViewerHDF5SetAIJNames()
418058bd781SVaclav Hapla @*/
419058bd781SVaclav Hapla PetscErrorCode  PetscViewerHDF5GetAIJNames(PetscViewer viewer, const char *iname[], const char *jname[], const char *aname[], const char *cname[])
420058bd781SVaclav Hapla {
421058bd781SVaclav Hapla   PetscErrorCode ierr;
422058bd781SVaclav Hapla 
423058bd781SVaclav Hapla   PetscFunctionBegin;
424058bd781SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
425058bd781SVaclav Hapla   PetscValidPointer(iname,2);
426058bd781SVaclav Hapla   PetscValidPointer(jname,3);
427058bd781SVaclav Hapla   PetscValidPointer(aname,4);
428058bd781SVaclav Hapla   PetscValidPointer(cname,5);
429058bd781SVaclav Hapla   ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetAIJNames_C",(PetscViewer,const char*[],const char*[],const char*[],const char*[]),(viewer,iname,jname,aname,cname));CHKERRQ(ierr);
430058bd781SVaclav Hapla   PetscFunctionReturn(0);
431058bd781SVaclav Hapla }
432058bd781SVaclav Hapla 
433b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
434b723ab35SVaclav Hapla {
435b723ab35SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
436b723ab35SVaclav Hapla   PetscErrorCode   ierr;
437b723ab35SVaclav Hapla 
438b723ab35SVaclav Hapla   PetscFunctionBegin;
439b723ab35SVaclav Hapla   if (!hdf5->mataij_names_set) {
440cbb4c999SVaclav Hapla     if (viewer->format == PETSC_VIEWER_HDF5_MAT) {
441b723ab35SVaclav Hapla       ierr = PetscViewerHDF5SetAIJNames_HDF5(viewer,"jc","ir","data","MATLAB_sparse");CHKERRQ(ierr);
442cbb4c999SVaclav Hapla     } else {
443cbb4c999SVaclav Hapla       ierr = PetscViewerHDF5SetAIJNames_HDF5(viewer,"i","j","a","ncols");CHKERRQ(ierr);
444cbb4c999SVaclav Hapla     }
445b723ab35SVaclav Hapla   }
446b723ab35SVaclav Hapla   PetscFunctionReturn(0);
447b723ab35SVaclav Hapla }
448b723ab35SVaclav Hapla 
4498556b5ebSBarry Smith /*MC
4508556b5ebSBarry Smith    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
4518556b5ebSBarry Smith 
4528556b5ebSBarry Smith 
4538556b5ebSBarry Smith .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
4548556b5ebSBarry Smith            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
4558556b5ebSBarry Smith            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
4568556b5ebSBarry Smith            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
4578556b5ebSBarry Smith 
4581b266c99SBarry Smith   Level: beginner
4598556b5ebSBarry Smith M*/
460d1232d7fSToby Isaac 
4618cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
4625c6c1daeSBarry Smith {
4635c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
4645c6c1daeSBarry Smith   PetscErrorCode   ierr;
4655c6c1daeSBarry Smith 
4665c6c1daeSBarry Smith   PetscFunctionBegin;
467b00a9115SJed Brown   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
4685c6c1daeSBarry Smith 
4695c6c1daeSBarry Smith   v->data                = (void*) hdf5;
4705c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
47182ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
472b723ab35SVaclav Hapla   v->ops->setup          = PetscViewerSetUp_HDF5;
4735c6c1daeSBarry Smith   v->ops->flush          = 0;
4745c6c1daeSBarry Smith   hdf5->btype            = (PetscFileMode) -1;
4755c6c1daeSBarry Smith   hdf5->filename         = 0;
4765c6c1daeSBarry Smith   hdf5->timestep         = -1;
4770298fd71SBarry Smith   hdf5->groups           = NULL;
4785c6c1daeSBarry Smith 
479b723ab35SVaclav Hapla   hdf5->mataij_iname     = NULL;
480b723ab35SVaclav Hapla   hdf5->mataij_jname     = NULL;
481b723ab35SVaclav Hapla   hdf5->mataij_aname     = NULL;
482b723ab35SVaclav Hapla   hdf5->mataij_cname     = NULL;
483b723ab35SVaclav Hapla   hdf5->mataij_names_set = PETSC_FALSE;
484058bd781SVaclav Hapla 
4850b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
486d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr);
4870b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
488*0b62783dSVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr);
48982ea9b62SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
4909a0502c6SHåkon Strandenes   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr);
491058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetAIJNames_C",PetscViewerHDF5SetAIJNames_HDF5);CHKERRQ(ierr);
492058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetAIJNames_C",PetscViewerHDF5GetAIJNames_HDF5);CHKERRQ(ierr);
4935c6c1daeSBarry Smith   PetscFunctionReturn(0);
4945c6c1daeSBarry Smith }
4955c6c1daeSBarry Smith 
4965c6c1daeSBarry Smith /*@C
4975c6c1daeSBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
4985c6c1daeSBarry Smith 
4995c6c1daeSBarry Smith    Collective on MPI_Comm
5005c6c1daeSBarry Smith 
5015c6c1daeSBarry Smith    Input Parameters:
5025c6c1daeSBarry Smith +  comm - MPI communicator
5035c6c1daeSBarry Smith .  name - name of file
5045c6c1daeSBarry Smith -  type - type of file
5055c6c1daeSBarry Smith $    FILE_MODE_WRITE - create new file for binary output
5065c6c1daeSBarry Smith $    FILE_MODE_READ - open existing file for binary input
5075c6c1daeSBarry Smith $    FILE_MODE_APPEND - open existing file for binary output
5085c6c1daeSBarry Smith 
5095c6c1daeSBarry Smith    Output Parameter:
5105c6c1daeSBarry Smith .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
5115c6c1daeSBarry Smith 
51282ea9b62SBarry Smith   Options Database:
51382ea9b62SBarry 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
5149a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
51582ea9b62SBarry Smith 
5165c6c1daeSBarry Smith    Level: beginner
5175c6c1daeSBarry Smith 
5185c6c1daeSBarry Smith    Note:
5195c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
5205c6c1daeSBarry Smith 
5215c6c1daeSBarry Smith    Concepts: HDF5 files
5225c6c1daeSBarry Smith    Concepts: PetscViewerHDF5^creating
5235c6c1daeSBarry Smith 
5246a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
5259a0502c6SHåkon Strandenes           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
526a56f64adSBarry Smith           MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
5275c6c1daeSBarry Smith @*/
5285c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
5295c6c1daeSBarry Smith {
5305c6c1daeSBarry Smith   PetscErrorCode ierr;
5315c6c1daeSBarry Smith 
5325c6c1daeSBarry Smith   PetscFunctionBegin;
5335c6c1daeSBarry Smith   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
5345c6c1daeSBarry Smith   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
5355c6c1daeSBarry Smith   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
5365c6c1daeSBarry Smith   ierr = PetscViewerFileSetName(*hdf5v, name);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 
576874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
5775c6c1daeSBarry Smith @*/
5785c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
5795c6c1daeSBarry Smith {
5805c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5815c6c1daeSBarry Smith   GroupList        *groupNode;
5825c6c1daeSBarry Smith   PetscErrorCode   ierr;
5835c6c1daeSBarry Smith 
5845c6c1daeSBarry Smith   PetscFunctionBegin;
5855c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5865c6c1daeSBarry Smith   PetscValidCharPointer(name,2);
58795dccacaSBarry Smith   ierr = PetscNew(&groupNode);CHKERRQ(ierr);
5885c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
589a297a907SKarl Rupp 
5905c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
5915c6c1daeSBarry Smith   hdf5->groups    = groupNode;
5925c6c1daeSBarry Smith   PetscFunctionReturn(0);
5935c6c1daeSBarry Smith }
5945c6c1daeSBarry Smith 
5953ef9c667SSatish Balay /*@
5965c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
5975c6c1daeSBarry Smith 
5985c6c1daeSBarry Smith   Not collective
5995c6c1daeSBarry Smith 
6005c6c1daeSBarry Smith   Input Parameter:
6015c6c1daeSBarry Smith . viewer - the PetscViewer
6025c6c1daeSBarry Smith 
6035c6c1daeSBarry Smith   Level: intermediate
6045c6c1daeSBarry Smith 
605874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
6065c6c1daeSBarry Smith @*/
6075c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
6085c6c1daeSBarry Smith {
6095c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
6105c6c1daeSBarry Smith   GroupList        *groupNode;
6115c6c1daeSBarry Smith   PetscErrorCode   ierr;
6125c6c1daeSBarry Smith 
6135c6c1daeSBarry Smith   PetscFunctionBegin;
6145c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
61582f516ccSBarry Smith   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
6165c6c1daeSBarry Smith   groupNode    = hdf5->groups;
6175c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
6185c6c1daeSBarry Smith   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
6195c6c1daeSBarry Smith   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
6205c6c1daeSBarry Smith   PetscFunctionReturn(0);
6215c6c1daeSBarry Smith }
6225c6c1daeSBarry Smith 
6235c6c1daeSBarry Smith /*@C
624874270d9SVaclav Hapla   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup().
625874270d9SVaclav Hapla   If none has been assigned, returns NULL.
6265c6c1daeSBarry Smith 
6275c6c1daeSBarry Smith   Not collective
6285c6c1daeSBarry Smith 
6295c6c1daeSBarry Smith   Input Parameter:
6305c6c1daeSBarry Smith . viewer - the PetscViewer
6315c6c1daeSBarry Smith 
6325c6c1daeSBarry Smith   Output Parameter:
6335c6c1daeSBarry Smith . name - The group name
6345c6c1daeSBarry Smith 
6355c6c1daeSBarry Smith   Level: intermediate
6365c6c1daeSBarry Smith 
637874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
6385c6c1daeSBarry Smith @*/
6395c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
6405c6c1daeSBarry Smith {
6415c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
6425c6c1daeSBarry Smith 
6435c6c1daeSBarry Smith   PetscFunctionBegin;
6445c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
645c959eef4SJed Brown   PetscValidPointer(name,2);
646a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
6470298fd71SBarry Smith   else *name = NULL;
6485c6c1daeSBarry Smith   PetscFunctionReturn(0);
6495c6c1daeSBarry Smith }
6505c6c1daeSBarry Smith 
6518c557b5aSVaclav Hapla /*@
652874270d9SVaclav Hapla   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(),
653874270d9SVaclav Hapla   and return this group's ID and file ID.
654874270d9SVaclav Hapla   If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID.
655874270d9SVaclav Hapla 
656874270d9SVaclav Hapla   Not collective
657874270d9SVaclav Hapla 
658874270d9SVaclav Hapla   Input Parameter:
659874270d9SVaclav Hapla . viewer - the PetscViewer
660874270d9SVaclav Hapla 
661874270d9SVaclav Hapla   Output Parameter:
662874270d9SVaclav Hapla + fileId - The HDF5 file ID
663874270d9SVaclav Hapla - groupId - The HDF5 group ID
664874270d9SVaclav Hapla 
665874270d9SVaclav Hapla   Level: intermediate
666874270d9SVaclav Hapla 
667874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
668874270d9SVaclav Hapla @*/
66954dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
67054dbf706SVaclav Hapla {
67190e03537SVaclav Hapla   hid_t          file_id;
67290e03537SVaclav Hapla   H5O_type_t     type;
67354dbf706SVaclav Hapla   const char     *groupName = NULL;
67454dbf706SVaclav Hapla   PetscErrorCode ierr;
67554dbf706SVaclav Hapla 
67654dbf706SVaclav Hapla   PetscFunctionBegin;
67754dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr);
67854dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr);
67990e03537SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, PETSC_TRUE, NULL, &type);CHKERRQ(ierr);
68090e03537SVaclav 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);
68190e03537SVaclav Hapla   PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT));
68254dbf706SVaclav Hapla   *fileId  = file_id;
68354dbf706SVaclav Hapla   PetscFunctionReturn(0);
68454dbf706SVaclav Hapla }
68554dbf706SVaclav Hapla 
6863ef9c667SSatish Balay /*@
6875c6c1daeSBarry Smith   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
6885c6c1daeSBarry Smith 
6895c6c1daeSBarry Smith   Not collective
6905c6c1daeSBarry Smith 
6915c6c1daeSBarry Smith   Input Parameter:
6925c6c1daeSBarry Smith . viewer - the PetscViewer
6935c6c1daeSBarry Smith 
6945c6c1daeSBarry Smith   Level: intermediate
6955c6c1daeSBarry Smith 
6965c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
6975c6c1daeSBarry Smith @*/
6985c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
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;
7055c6c1daeSBarry Smith   PetscFunctionReturn(0);
7065c6c1daeSBarry Smith }
7075c6c1daeSBarry Smith 
7083ef9c667SSatish Balay /*@
7095c6c1daeSBarry Smith   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
7105c6c1daeSBarry Smith   of -1 disables blocking with timesteps.
7115c6c1daeSBarry Smith 
7125c6c1daeSBarry Smith   Not collective
7135c6c1daeSBarry Smith 
7145c6c1daeSBarry Smith   Input Parameters:
7155c6c1daeSBarry Smith + viewer - the PetscViewer
7165c6c1daeSBarry Smith - timestep - The timestep number
7175c6c1daeSBarry Smith 
7185c6c1daeSBarry Smith   Level: intermediate
7195c6c1daeSBarry Smith 
7205c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
7215c6c1daeSBarry Smith @*/
7225c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
7235c6c1daeSBarry Smith {
7245c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7255c6c1daeSBarry Smith 
7265c6c1daeSBarry Smith   PetscFunctionBegin;
7275c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7285c6c1daeSBarry Smith   hdf5->timestep = timestep;
7295c6c1daeSBarry Smith   PetscFunctionReturn(0);
7305c6c1daeSBarry Smith }
7315c6c1daeSBarry Smith 
7323ef9c667SSatish Balay /*@
7335c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
7345c6c1daeSBarry Smith 
7355c6c1daeSBarry Smith   Not collective
7365c6c1daeSBarry Smith 
7375c6c1daeSBarry Smith   Input Parameter:
7385c6c1daeSBarry Smith . viewer - the PetscViewer
7395c6c1daeSBarry Smith 
7405c6c1daeSBarry Smith   Output Parameter:
7415c6c1daeSBarry Smith . timestep - The timestep number
7425c6c1daeSBarry Smith 
7435c6c1daeSBarry Smith   Level: intermediate
7445c6c1daeSBarry Smith 
7455c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
7465c6c1daeSBarry Smith @*/
7475c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
7485c6c1daeSBarry Smith {
7495c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7505c6c1daeSBarry Smith 
7515c6c1daeSBarry Smith   PetscFunctionBegin;
7525c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7535c6c1daeSBarry Smith   PetscValidPointer(timestep,2);
7545c6c1daeSBarry Smith   *timestep = hdf5->timestep;
7555c6c1daeSBarry Smith   PetscFunctionReturn(0);
7565c6c1daeSBarry Smith }
7575c6c1daeSBarry Smith 
75836bce6e8SMatthew G. Knepley /*@C
75936bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
76036bce6e8SMatthew G. Knepley 
76136bce6e8SMatthew G. Knepley   Not collective
76236bce6e8SMatthew G. Knepley 
76336bce6e8SMatthew G. Knepley   Input Parameter:
76436bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
76536bce6e8SMatthew G. Knepley 
76636bce6e8SMatthew G. Knepley   Output Parameter:
76736bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
76836bce6e8SMatthew G. Knepley 
76936bce6e8SMatthew G. Knepley   Level: advanced
77036bce6e8SMatthew G. Knepley 
77136bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
77236bce6e8SMatthew G. Knepley @*/
77336bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
77436bce6e8SMatthew G. Knepley {
77536bce6e8SMatthew G. Knepley   PetscFunctionBegin;
77636bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
77736bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
77836bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_LLONG;
77936bce6e8SMatthew G. Knepley #else
78036bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_INT;
78136bce6e8SMatthew G. Knepley #endif
78236bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
78336bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
78436bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
78536bce6e8SMatthew G. Knepley   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_DOUBLE;
786de3ba262SVaclav Hapla   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_INT;
78736bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
78836bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
78936bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
7907e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
79136bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
79236bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
79336bce6e8SMatthew G. Knepley }
79436bce6e8SMatthew G. Knepley 
79536bce6e8SMatthew G. Knepley /*@C
79636bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
79736bce6e8SMatthew G. Knepley 
79836bce6e8SMatthew G. Knepley   Not collective
79936bce6e8SMatthew G. Knepley 
80036bce6e8SMatthew G. Knepley   Input Parameter:
80136bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
80236bce6e8SMatthew G. Knepley 
80336bce6e8SMatthew G. Knepley   Output Parameter:
80436bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
80536bce6e8SMatthew G. Knepley 
80636bce6e8SMatthew G. Knepley   Level: advanced
80736bce6e8SMatthew G. Knepley 
80836bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
80936bce6e8SMatthew G. Knepley @*/
81036bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
81136bce6e8SMatthew G. Knepley {
81236bce6e8SMatthew G. Knepley   PetscFunctionBegin;
81336bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
81436bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
81536bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
81636bce6e8SMatthew G. Knepley #else
81736bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
81836bce6e8SMatthew G. Knepley #endif
81936bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
82036bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
82136bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
82236bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
82336bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
82436bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
8257e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
82636bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
82736bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
82836bce6e8SMatthew G. Knepley }
82936bce6e8SMatthew G. Knepley 
830df863907SAlex Fikl /*@C
831b8ee85a1SVaclav Hapla  PetscViewerHDF5WriteAttribute - Write an attribute
83236bce6e8SMatthew G. Knepley 
83336bce6e8SMatthew G. Knepley   Input Parameters:
83436bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer
83557d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
83636bce6e8SMatthew G. Knepley . name   - The attribute name
83736bce6e8SMatthew G. Knepley . datatype - The attribute type
83836bce6e8SMatthew G. Knepley - value    - The attribute value
83936bce6e8SMatthew G. Knepley 
84036bce6e8SMatthew G. Knepley   Level: advanced
84136bce6e8SMatthew G. Knepley 
842577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
84336bce6e8SMatthew G. Knepley @*/
84457d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value)
84536bce6e8SMatthew G. Knepley {
84657d22f7dSVaclav Hapla   char           *parent;
84760568592SMatthew G. Knepley   hid_t          h5, dataspace, obj, attribute, dtype;
848080f144cSVaclav Hapla   PetscBool      has;
84936bce6e8SMatthew G. Knepley   PetscErrorCode ierr;
85036bce6e8SMatthew G. Knepley 
85136bce6e8SMatthew G. Knepley   PetscFunctionBegin;
8525cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
85357d22f7dSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset, 2);
854c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
855b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
85657d22f7dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
857bb286ee1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr);
858b67695ecSVaclav Hapla   ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);
85936bce6e8SMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
8607e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
8617e97c476SMatthew G. Knepley     size_t len;
8627e97c476SMatthew G. Knepley     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
863729ab6d9SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
8647e97c476SMatthew G. Knepley   }
86536bce6e8SMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
866729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
86760568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
868080f144cSVaclav Hapla   if (has) {
869080f144cSVaclav Hapla     PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
870080f144cSVaclav Hapla   } else {
87160568592SMatthew G. Knepley     PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
872080f144cSVaclav Hapla   }
873729ab6d9SBarry Smith   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
874729ab6d9SBarry Smith   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
875729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
87660568592SMatthew G. Knepley   PetscStackCallHDF5(H5Oclose,(obj));
877729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(dataspace));
87857d22f7dSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
87936bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
88036bce6e8SMatthew G. Knepley }
88136bce6e8SMatthew G. Knepley 
882df863907SAlex Fikl /*@C
883577e0e04SVaclav Hapla  PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name
884577e0e04SVaclav Hapla 
885577e0e04SVaclav Hapla   Input Parameters:
886577e0e04SVaclav Hapla + viewer   - The HDF5 viewer
887577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
888577e0e04SVaclav Hapla . name     - The attribute name
889577e0e04SVaclav Hapla . datatype - The attribute type
890577e0e04SVaclav Hapla - value    - The attribute value
891577e0e04SVaclav Hapla 
892577e0e04SVaclav Hapla   Notes:
893577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
894577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
895577e0e04SVaclav Hapla 
896577e0e04SVaclav Hapla   Level: advanced
897577e0e04SVaclav Hapla 
898577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
899577e0e04SVaclav Hapla @*/
900577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
901577e0e04SVaclav Hapla {
902577e0e04SVaclav Hapla   PetscErrorCode ierr;
903577e0e04SVaclav Hapla 
904577e0e04SVaclav Hapla   PetscFunctionBegin;
905577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
906577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
907577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
908b7e81f0fSVaclav Hapla   PetscValidPointer(value,5);
909577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
910577e0e04SVaclav Hapla   ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
911577e0e04SVaclav Hapla   PetscFunctionReturn(0);
912577e0e04SVaclav Hapla }
913577e0e04SVaclav Hapla 
914577e0e04SVaclav Hapla /*@C
915b8ee85a1SVaclav Hapla  PetscViewerHDF5ReadAttribute - Read an attribute
916ad0c61feSMatthew G. Knepley 
917ad0c61feSMatthew G. Knepley   Input Parameters:
918ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer
91957d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
920ad0c61feSMatthew G. Knepley . name   - The attribute name
921ad0c61feSMatthew G. Knepley - datatype - The attribute type
922ad0c61feSMatthew G. Knepley 
923ad0c61feSMatthew G. Knepley   Output Parameter:
924ad0c61feSMatthew G. Knepley . value    - The attribute value
925ad0c61feSMatthew G. Knepley 
926ad0c61feSMatthew G. Knepley   Level: advanced
927ad0c61feSMatthew G. Knepley 
928577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
929ad0c61feSMatthew G. Knepley @*/
93057d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value)
931ad0c61feSMatthew G. Knepley {
93257d22f7dSVaclav Hapla   char           *parent;
933f0b58479SMatthew G. Knepley   hid_t          h5, obj, attribute, atype, dtype;
934969748fdSVaclav Hapla   PetscBool      has;
935ad0c61feSMatthew G. Knepley   PetscErrorCode ierr;
936ad0c61feSMatthew G. Knepley 
937ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
9385cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
93957d22f7dSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset, 2);
940c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
941b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
94257d22f7dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
943969748fdSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr);
944969748fdSVaclav Hapla   if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);}
945969748fdSVaclav Hapla   if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name);
946ad0c61feSMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
947ad0c61feSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
94860568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
94960568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
950f0b58479SMatthew G. Knepley   if (datatype == PETSC_STRING) {
951f0b58479SMatthew G. Knepley     size_t len;
952e90c5075SMarek Pecha     PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
95315b861d2SVaclav Hapla     PetscStackCall("H5Tget_size",len = H5Tget_size(atype));
954f0b58479SMatthew G. Knepley     PetscStackCallHDF5(H5Tclose,(atype));
955f0b58479SMatthew G. Knepley     ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr);
956f0b58479SMatthew G. Knepley   }
95770efba33SBarry Smith   PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
958729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
959e90c5075SMarek Pecha   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
960e90c5075SMarek Pecha   PetscStackCallHDF5(H5Oclose,(obj));
96157d22f7dSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
962ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
963ad0c61feSMatthew G. Knepley }
964ad0c61feSMatthew G. Knepley 
965577e0e04SVaclav Hapla /*@C
966577e0e04SVaclav Hapla  PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name
967577e0e04SVaclav Hapla 
968577e0e04SVaclav Hapla   Input Parameters:
969577e0e04SVaclav Hapla + viewer   - The HDF5 viewer
970577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
971577e0e04SVaclav Hapla . name     - The attribute name
972577e0e04SVaclav Hapla - datatype - The attribute type
973577e0e04SVaclav Hapla 
974577e0e04SVaclav Hapla   Output Parameter:
975577e0e04SVaclav Hapla . value    - The attribute value
976577e0e04SVaclav Hapla 
977577e0e04SVaclav Hapla   Notes:
978577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
979577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
980577e0e04SVaclav Hapla 
981577e0e04SVaclav Hapla   Level: advanced
982577e0e04SVaclav Hapla 
983577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
984577e0e04SVaclav Hapla @*/
985577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value)
986577e0e04SVaclav Hapla {
987577e0e04SVaclav Hapla   PetscErrorCode ierr;
988577e0e04SVaclav Hapla 
989577e0e04SVaclav Hapla   PetscFunctionBegin;
990577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
991577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
992577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
993b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
994577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
995577e0e04SVaclav Hapla   ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
996577e0e04SVaclav Hapla   PetscFunctionReturn(0);
997577e0e04SVaclav Hapla }
998577e0e04SVaclav Hapla 
999a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
1000a75c3fd4SVaclav Hapla {
1001a75c3fd4SVaclav Hapla   htri_t exists;
1002a75c3fd4SVaclav Hapla   hid_t group;
1003a75c3fd4SVaclav Hapla 
1004a75c3fd4SVaclav Hapla   PetscFunctionBegin;
1005a75c3fd4SVaclav Hapla   PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT));
1006a75c3fd4SVaclav Hapla   if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT));
1007a75c3fd4SVaclav Hapla   if (!exists && createGroup) {
1008a75c3fd4SVaclav Hapla     PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1009a75c3fd4SVaclav Hapla     PetscStackCallHDF5(H5Gclose,(group));
1010a75c3fd4SVaclav Hapla     exists = PETSC_TRUE;
1011a75c3fd4SVaclav Hapla   }
1012a75c3fd4SVaclav Hapla   *exists_ = (PetscBool) exists;
1013a75c3fd4SVaclav Hapla   PetscFunctionReturn(0);
1014a75c3fd4SVaclav Hapla }
1015a75c3fd4SVaclav Hapla 
1016bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
10175cdeb108SMatthew G. Knepley {
101890e03537SVaclav Hapla   const char     rootGroupName[] = "/";
10195cdeb108SMatthew G. Knepley   hid_t          h5;
1020e5bf9ebcSVaclav Hapla   PetscBool      exists=PETSC_FALSE;
102115b861d2SVaclav Hapla   PetscInt       i;
102215b861d2SVaclav Hapla   int            n;
102385688ae2SVaclav Hapla   char           **hierarchy;
102485688ae2SVaclav Hapla   char           buf[PETSC_MAX_PATH_LEN]="";
10255cdeb108SMatthew G. Knepley   PetscErrorCode ierr;
10265cdeb108SMatthew G. Knepley 
10275cdeb108SMatthew G. Knepley   PetscFunctionBegin;
10285cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
102990e03537SVaclav Hapla   if (name) PetscValidCharPointer(name, 2);
103090e03537SVaclav Hapla   else name = rootGroupName;
1031ccd66a83SVaclav Hapla   if (has) {
103256cc0592SVaclav Hapla     PetscValidIntPointer(has, 3);
10335cdeb108SMatthew G. Knepley     *has = PETSC_FALSE;
1034ccd66a83SVaclav Hapla   }
1035ccd66a83SVaclav Hapla   if (otype) {
1036ccd66a83SVaclav Hapla     PetscValidIntPointer(otype, 4);
103756cc0592SVaclav Hapla     *otype = H5O_TYPE_UNKNOWN;
1038ccd66a83SVaclav Hapla   }
1039c0cd0301SJed Brown   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
104085688ae2SVaclav Hapla 
104185688ae2SVaclav Hapla   /*
104285688ae2SVaclav Hapla      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
104385688ae2SVaclav Hapla      Hence, each of them needs to be tested separately:
104485688ae2SVaclav Hapla      1) whether it's a valid link
104585688ae2SVaclav Hapla      2) whether this link resolves to an object
104685688ae2SVaclav Hapla      See H5Oexists_by_name() documentation.
104785688ae2SVaclav Hapla   */
104885688ae2SVaclav Hapla   ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr);
104985688ae2SVaclav Hapla   if (!n) {
105085688ae2SVaclav Hapla     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1051ccd66a83SVaclav Hapla     if (has)   *has   = PETSC_TRUE;
1052ccd66a83SVaclav Hapla     if (otype) *otype = H5O_TYPE_GROUP;
105385688ae2SVaclav Hapla     ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
105485688ae2SVaclav Hapla     PetscFunctionReturn(0);
105585688ae2SVaclav Hapla   }
105685688ae2SVaclav Hapla   for (i=0; i<n; i++) {
105785688ae2SVaclav Hapla     ierr = PetscStrcat(buf,"/");CHKERRQ(ierr);
105885688ae2SVaclav Hapla     ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr);
1059a75c3fd4SVaclav Hapla     ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr);
1060a75c3fd4SVaclav Hapla     if (!exists) break;
106185688ae2SVaclav Hapla   }
106285688ae2SVaclav Hapla   ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
106385688ae2SVaclav Hapla 
106485688ae2SVaclav Hapla   /* If the object exists, get its type */
1065ccd66a83SVaclav Hapla   if (exists && otype) {
10665cdeb108SMatthew G. Knepley     H5O_info_t info;
10675cdeb108SMatthew G. Knepley 
106876276748SVaclav Hapla     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
106904633f7fSVaclav Hapla     PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT));
107056cc0592SVaclav Hapla     *otype = info.type;
10715cdeb108SMatthew G. Knepley   }
1072ccd66a83SVaclav Hapla   if (has) *has = exists;
10735cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
10745cdeb108SMatthew G. Knepley }
10755cdeb108SMatthew G. Knepley 
1076ecce1506SVaclav Hapla /*@
10770a9f38efSVaclav Hapla  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
10780a9f38efSVaclav Hapla 
10790a9f38efSVaclav Hapla   Input Parameters:
10800a9f38efSVaclav Hapla . viewer - The HDF5 viewer
10810a9f38efSVaclav Hapla 
10820a9f38efSVaclav Hapla   Output Parameter:
10830a9f38efSVaclav Hapla . has    - Flag for group existence
10840a9f38efSVaclav Hapla 
10850a9f38efSVaclav Hapla   Notes:
10860a9f38efSVaclav Hapla   If the path exists but is not a group, this returns PETSC_FALSE as well.
10870a9f38efSVaclav Hapla 
10880a9f38efSVaclav Hapla   Level: advanced
10890a9f38efSVaclav Hapla 
10900a9f38efSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup()
10910a9f38efSVaclav Hapla @*/
10920a9f38efSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has)
10930a9f38efSVaclav Hapla {
10940a9f38efSVaclav Hapla   H5O_type_t type;
10950a9f38efSVaclav Hapla   const char *name;
10960a9f38efSVaclav Hapla   PetscErrorCode ierr;
10970a9f38efSVaclav Hapla 
10980a9f38efSVaclav Hapla   PetscFunctionBegin;
10990a9f38efSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
11000a9f38efSVaclav Hapla   PetscValidIntPointer(has,2);
11010a9f38efSVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr);
11020a9f38efSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr);
11030a9f38efSVaclav Hapla   *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE;
11040a9f38efSVaclav Hapla   PetscFunctionReturn(0);
11050a9f38efSVaclav Hapla }
11060a9f38efSVaclav Hapla 
11070a9f38efSVaclav Hapla /*@
1108e3f143f7SVaclav Hapla  PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1109ecce1506SVaclav Hapla 
1110ecce1506SVaclav Hapla   Input Parameters:
1111ecce1506SVaclav Hapla + viewer - The HDF5 viewer
1112ecce1506SVaclav Hapla - obj    - The named object
1113ecce1506SVaclav Hapla 
1114ecce1506SVaclav Hapla   Output Parameter:
1115ecce1506SVaclav Hapla . has    - Flag for dataset existence; PETSC_FALSE for unnamed object
1116ecce1506SVaclav Hapla 
1117e3f143f7SVaclav Hapla   Notes:
1118e3f143f7SVaclav Hapla   If the path exists but is not a dataset, this returns PETSC_FALSE as well.
1119e3f143f7SVaclav Hapla 
1120ecce1506SVaclav Hapla   Level: advanced
1121ecce1506SVaclav Hapla 
1122e3f143f7SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1123ecce1506SVaclav Hapla @*/
1124ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1125ecce1506SVaclav Hapla {
112656cc0592SVaclav Hapla   H5O_type_t type;
11276c132bc1SVaclav Hapla   char *path;
1128ecce1506SVaclav Hapla   PetscErrorCode ierr;
1129ecce1506SVaclav Hapla 
1130ecce1506SVaclav Hapla   PetscFunctionBegin;
1131c57a1a9aSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1132c57a1a9aSVaclav Hapla   PetscValidHeader(obj,2);
1133c57a1a9aSVaclav Hapla   PetscValidIntPointer(has,3);
1134ecce1506SVaclav Hapla   *has = PETSC_FALSE;
1135e3f143f7SVaclav Hapla   if (!obj->name) PetscFunctionReturn(0);
11366c132bc1SVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr);
11376c132bc1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr);
113856cc0592SVaclav Hapla   *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE;
11396c132bc1SVaclav Hapla   ierr = PetscFree(path);CHKERRQ(ierr);
1140ecce1506SVaclav Hapla   PetscFunctionReturn(0);
1141ecce1506SVaclav Hapla }
1142ecce1506SVaclav Hapla 
1143df863907SAlex Fikl /*@C
1144ece83b6cSVaclav Hapla  PetscViewerHDF5HasAttribute - Check whether an attribute exists
1145e7bdbf8eSMatthew G. Knepley 
1146e7bdbf8eSMatthew G. Knepley   Input Parameters:
1147e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer
114810e69e8fSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
1149e7bdbf8eSMatthew G. Knepley - name   - The attribute name
1150e7bdbf8eSMatthew G. Knepley 
1151e7bdbf8eSMatthew G. Knepley   Output Parameter:
1152e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
1153e7bdbf8eSMatthew G. Knepley 
1154e7bdbf8eSMatthew G. Knepley   Level: advanced
1155e7bdbf8eSMatthew G. Knepley 
1156577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1157e7bdbf8eSMatthew G. Knepley @*/
115810e69e8fSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has)
1159e7bdbf8eSMatthew G. Knepley {
116010e69e8fSVaclav Hapla   char           *parent;
1161e7bdbf8eSMatthew G. Knepley   PetscErrorCode ierr;
1162e7bdbf8eSMatthew G. Knepley 
1163e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
11645cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
116510e69e8fSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset,2);
1166c57a1a9aSVaclav Hapla   PetscValidCharPointer(name,3);
1167c57a1a9aSVaclav Hapla   PetscValidIntPointer(has,4);
116810e69e8fSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
1169bb286ee1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr);
117010e69e8fSVaclav Hapla   if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);}
117110e69e8fSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
117206db490cSVaclav Hapla   PetscFunctionReturn(0);
117306db490cSVaclav Hapla }
117406db490cSVaclav Hapla 
1175577e0e04SVaclav Hapla /*@C
1176577e0e04SVaclav Hapla  PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name
1177577e0e04SVaclav Hapla 
1178577e0e04SVaclav Hapla   Input Parameters:
1179577e0e04SVaclav Hapla + viewer - The HDF5 viewer
1180577e0e04SVaclav Hapla . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1181577e0e04SVaclav Hapla - name   - The attribute name
1182577e0e04SVaclav Hapla 
1183577e0e04SVaclav Hapla   Output Parameter:
1184577e0e04SVaclav Hapla . has    - Flag for attribute existence
1185577e0e04SVaclav Hapla 
1186577e0e04SVaclav Hapla   Notes:
1187577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1188577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
1189577e0e04SVaclav Hapla 
1190577e0e04SVaclav Hapla   Level: advanced
1191577e0e04SVaclav Hapla 
1192577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1193577e0e04SVaclav Hapla @*/
1194577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1195577e0e04SVaclav Hapla {
1196577e0e04SVaclav Hapla   PetscErrorCode ierr;
1197577e0e04SVaclav Hapla 
1198577e0e04SVaclav Hapla   PetscFunctionBegin;
1199577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1200577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
1201577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
1202577e0e04SVaclav Hapla   PetscValidIntPointer(has,4);
1203577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
1204577e0e04SVaclav Hapla   ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr);
1205577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1206577e0e04SVaclav Hapla }
1207577e0e04SVaclav Hapla 
120806db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
120906db490cSVaclav Hapla {
121006db490cSVaclav Hapla   hid_t          h5;
121106db490cSVaclav Hapla   htri_t         hhas;
121206db490cSVaclav Hapla   PetscErrorCode ierr;
121306db490cSVaclav Hapla 
121406db490cSVaclav Hapla   PetscFunctionBegin;
121506db490cSVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
12162f936e54SVaclav Hapla   PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT));
12175cdeb108SMatthew G. Knepley   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1218e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
1219e7bdbf8eSMatthew G. Knepley }
1220e7bdbf8eSMatthew G. Knepley 
1221eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadInitialize_Private(PetscViewer viewer, const char name[], HDF5ReadCtx *ctx)
1222a5e1feadSVaclav Hapla {
1223fbd37863SVaclav Hapla   HDF5ReadCtx    h=NULL;
1224a5e1feadSVaclav Hapla   PetscErrorCode ierr;
1225a5e1feadSVaclav Hapla 
1226a5e1feadSVaclav Hapla   PetscFunctionBegin;
122769a06e7bSVaclav Hapla   ierr = PetscNew(&h);CHKERRQ(ierr);
122869a06e7bSVaclav Hapla   ierr = PetscViewerHDF5OpenGroup(viewer, &h->file, &h->group);CHKERRQ(ierr);
122969a06e7bSVaclav Hapla   PetscStackCallHDF5Return(h->dataset,H5Dopen2,(h->group, name, H5P_DEFAULT));
123069a06e7bSVaclav Hapla   PetscStackCallHDF5Return(h->dataspace,H5Dget_space,(h->dataset));
12319568d583SVaclav Hapla   ierr = PetscViewerHDF5GetTimestep(viewer, &h->timestep);CHKERRQ(ierr);
1232adb363a2SVaclav Hapla   ierr = PetscViewerHDF5HasAttribute(viewer,name,"complex",&h->complexVal);CHKERRQ(ierr);
1233aa54fc5aSVaclav Hapla   if (h->complexVal) {ierr = PetscViewerHDF5ReadAttribute(viewer,name,"complex",PETSC_BOOL,&h->complexVal);CHKERRQ(ierr);}
123409dabeb0SVaclav Hapla   /* MATLAB stores column vectors horizontally */
1235adb363a2SVaclav Hapla   ierr = PetscViewerHDF5HasAttribute(viewer,name,"MATLAB_class",&h->horizontal);CHKERRQ(ierr);
1236b8ef406cSVaclav Hapla   /* Create property list for collective dataset read */
1237b8ef406cSVaclav Hapla   PetscStackCallHDF5Return(h->plist,H5Pcreate,(H5P_DATASET_XFER));
1238b8ef406cSVaclav Hapla #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
1239b8ef406cSVaclav Hapla   PetscStackCallHDF5(H5Pset_dxpl_mpio,(h->plist, H5FD_MPIO_COLLECTIVE));
1240b8ef406cSVaclav Hapla #endif
1241b8ef406cSVaclav Hapla   /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */
124269a06e7bSVaclav Hapla   *ctx = h;
124369a06e7bSVaclav Hapla   PetscFunctionReturn(0);
124469a06e7bSVaclav Hapla }
124569a06e7bSVaclav Hapla 
1246eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadFinalize_Private(PetscViewer viewer, HDF5ReadCtx *ctx)
124769a06e7bSVaclav Hapla {
124869a06e7bSVaclav Hapla   HDF5ReadCtx    h;
124969a06e7bSVaclav Hapla   PetscErrorCode ierr;
125069a06e7bSVaclav Hapla 
125169a06e7bSVaclav Hapla   PetscFunctionBegin;
125269a06e7bSVaclav Hapla   h = *ctx;
1253b8ef406cSVaclav Hapla   PetscStackCallHDF5(H5Pclose,(h->plist));
125490e03537SVaclav Hapla   PetscStackCallHDF5(H5Gclose,(h->group));
125569a06e7bSVaclav Hapla   PetscStackCallHDF5(H5Sclose,(h->dataspace));
125669a06e7bSVaclav Hapla   PetscStackCallHDF5(H5Dclose,(h->dataset));
125769a06e7bSVaclav Hapla   ierr = PetscFree(*ctx);CHKERRQ(ierr);
125869a06e7bSVaclav Hapla   PetscFunctionReturn(0);
125969a06e7bSVaclav Hapla }
126069a06e7bSVaclav Hapla 
1261eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadSizes_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout *map_)
126269a06e7bSVaclav Hapla {
126369a06e7bSVaclav Hapla   int            rdim, dim;
126469a06e7bSVaclav Hapla   hsize_t        dims[4];
126509dabeb0SVaclav Hapla   PetscInt       bsInd, lenInd, bs, len, N;
12668374c777SVaclav Hapla   PetscLayout    map;
12678374c777SVaclav Hapla   PetscErrorCode ierr;
126869a06e7bSVaclav Hapla 
126969a06e7bSVaclav Hapla   PetscFunctionBegin;
12708374c777SVaclav Hapla   if (!(*map_)) {
12718374c777SVaclav Hapla     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)viewer),map_);CHKERRQ(ierr);
12728374c777SVaclav Hapla   }
12738374c777SVaclav Hapla   map = *map_;
12749787f08bSVaclav Hapla   /* calculate expected number of dimensions */
1275a5e1feadSVaclav Hapla   dim = 0;
12769568d583SVaclav Hapla   if (ctx->timestep >= 0) ++dim;
1277a5e1feadSVaclav Hapla   ++dim; /* length in blocks */
12789568d583SVaclav Hapla   if (ctx->complexVal) ++dim;
12799787f08bSVaclav Hapla   /* get actual number of dimensions in dataset */
128069a06e7bSVaclav Hapla   PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(ctx->dataspace, dims, NULL));
12819787f08bSVaclav Hapla   /* calculate expected dimension indices */
12829787f08bSVaclav Hapla   lenInd = 0;
12839568d583SVaclav Hapla   if (ctx->timestep >= 0) ++lenInd;
12849787f08bSVaclav Hapla   bsInd = lenInd + 1;
12859568d583SVaclav Hapla   ctx->dim2 = PETSC_FALSE;
12869787f08bSVaclav Hapla   if (rdim == dim) {
128745e21b7fSVaclav Hapla     bs = 1; /* support vectors stored as 1D array */
12889787f08bSVaclav Hapla   } else if (rdim == dim+1) {
128945e21b7fSVaclav Hapla     bs = (PetscInt) dims[bsInd];
12909568d583SVaclav Hapla     if (bs == 1) ctx->dim2 = PETSC_TRUE; /* vector with blocksize of 1, still stored as 2D array */
1291a5e1feadSVaclav Hapla   } else {
1292b54f1188SVaclav Hapla     SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Number of dimensions %d not %d as expected", rdim, dim);
1293a5e1feadSVaclav Hapla   }
129409dabeb0SVaclav Hapla   len = dims[lenInd];
129509dabeb0SVaclav Hapla   if (ctx->horizontal) {
1296b54f1188SVaclav Hapla     if (len != 1) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Cannot have horizontal array with number of rows > 1. In case of MATLAB MAT-file, vectors must be saved as column vectors.");
129709dabeb0SVaclav Hapla     len = bs;
129809dabeb0SVaclav Hapla     bs = 1;
129909dabeb0SVaclav Hapla   }
130009dabeb0SVaclav Hapla   N = (PetscInt) len*bs;
13018374c777SVaclav Hapla 
13028374c777SVaclav Hapla   /* Set Vec sizes,blocksize,and type if not already set */
13038374c777SVaclav Hapla   if (map->bs < 0) {
130445e21b7fSVaclav Hapla     ierr = PetscLayoutSetBlockSize(map, bs);CHKERRQ(ierr);
130545e21b7fSVaclav Hapla   } else if (map->bs != bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Block size of array in file is %D, not %D as expected",bs,map->bs);
13068374c777SVaclav Hapla   if (map->N < 0) {
130745e21b7fSVaclav Hapla     ierr = PetscLayoutSetSize(map, N);CHKERRQ(ierr);
1308b54f1188SVaclav Hapla   } else if (map->N != N) SETERRQ2(PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED, "Global size of array in file is %D, not %D as expected",N,map->N);
130969a06e7bSVaclav Hapla   PetscFunctionReturn(0);
131069a06e7bSVaclav Hapla }
131169a06e7bSVaclav Hapla 
1312eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadSelectHyperslab_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout map, hid_t *memspace)
131316f6402dSVaclav Hapla {
131416f6402dSVaclav Hapla   hsize_t        count[4], offset[4];
131516f6402dSVaclav Hapla   int            dim;
131616f6402dSVaclav Hapla   PetscInt       bs, n, low;
131716f6402dSVaclav Hapla   PetscErrorCode ierr;
131816f6402dSVaclav Hapla 
131916f6402dSVaclav Hapla   PetscFunctionBegin;
132016f6402dSVaclav Hapla   /* Compute local size and ownership range */
132116f6402dSVaclav Hapla   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
132216f6402dSVaclav Hapla   ierr = PetscLayoutGetBlockSize(map, &bs);CHKERRQ(ierr);
132316f6402dSVaclav Hapla   ierr = PetscLayoutGetLocalSize(map, &n);CHKERRQ(ierr);
132416f6402dSVaclav Hapla   ierr = PetscLayoutGetRange(map, &low, NULL);CHKERRQ(ierr);
132516f6402dSVaclav Hapla 
132616f6402dSVaclav Hapla   /* Each process defines a dataset and reads it from the hyperslab in the file */
132716f6402dSVaclav Hapla   dim  = 0;
132816f6402dSVaclav Hapla   if (ctx->timestep >= 0) {
132916f6402dSVaclav Hapla     count[dim]  = 1;
133016f6402dSVaclav Hapla     offset[dim] = ctx->timestep;
133116f6402dSVaclav Hapla     ++dim;
133216f6402dSVaclav Hapla   }
133309dabeb0SVaclav Hapla   if (ctx->horizontal) {
133409dabeb0SVaclav Hapla     count[dim]  = 1;
133509dabeb0SVaclav Hapla     offset[dim] = 0;
133609dabeb0SVaclav Hapla     ++dim;
133709dabeb0SVaclav Hapla   }
133816f6402dSVaclav Hapla   {
133916f6402dSVaclav Hapla     ierr = PetscHDF5IntCast(n/bs, &count[dim]);CHKERRQ(ierr);
134016f6402dSVaclav Hapla     ierr = PetscHDF5IntCast(low/bs, &offset[dim]);CHKERRQ(ierr);
134116f6402dSVaclav Hapla     ++dim;
134216f6402dSVaclav Hapla   }
134316f6402dSVaclav Hapla   if (bs > 1 || ctx->dim2) {
134409dabeb0SVaclav Hapla     if (PetscUnlikely(ctx->horizontal)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "cannot have horizontal array with blocksize > 1");
134516f6402dSVaclav Hapla     count[dim]  = bs;
134616f6402dSVaclav Hapla     offset[dim] = 0;
134716f6402dSVaclav Hapla     ++dim;
134816f6402dSVaclav Hapla   }
134916f6402dSVaclav Hapla   if (ctx->complexVal) {
135016f6402dSVaclav Hapla     count[dim]  = 2;
135116f6402dSVaclav Hapla     offset[dim] = 0;
135216f6402dSVaclav Hapla     ++dim;
135316f6402dSVaclav Hapla   }
135416f6402dSVaclav Hapla   PetscStackCallHDF5Return(*memspace,H5Screate_simple,(dim, count, NULL));
135516f6402dSVaclav Hapla   PetscStackCallHDF5(H5Sselect_hyperslab,(ctx->dataspace, H5S_SELECT_SET, offset, NULL, count, NULL));
135616f6402dSVaclav Hapla   PetscFunctionReturn(0);
135716f6402dSVaclav Hapla }
135816f6402dSVaclav Hapla 
1359eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadArray_Private(PetscViewer viewer, HDF5ReadCtx h, hid_t datatype, hid_t memspace, void *arr)
1360ef2d82ceSVaclav Hapla {
1361ef2d82ceSVaclav Hapla   PetscFunctionBegin;
1362ef2d82ceSVaclav Hapla   PetscStackCallHDF5(H5Dread,(h->dataset, datatype, memspace, h->dataspace, h->plist, arr));
1363ef2d82ceSVaclav Hapla   PetscFunctionReturn(0);
1364ef2d82ceSVaclav Hapla }
1365ef2d82ceSVaclav Hapla 
1366dad982a8SVaclav Hapla PetscErrorCode PetscViewerHDF5Load(PetscViewer viewer, const char *name, PetscLayout map, hid_t datatype, void **newarr)
1367a25c73c6SVaclav Hapla {
1368fbd37863SVaclav Hapla   HDF5ReadCtx     h=NULL;
1369fbd37863SVaclav Hapla   hid_t           memspace=0;
1370a25c73c6SVaclav Hapla   size_t          unitsize;
1371a25c73c6SVaclav Hapla   void            *arr;
1372a25c73c6SVaclav Hapla   PetscErrorCode  ierr;
1373a25c73c6SVaclav Hapla 
1374a25c73c6SVaclav Hapla   PetscFunctionBegin;
1375eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr);
13765a89fdf4SVaclav Hapla #if defined(PETSC_USE_COMPLEX)
13775a89fdf4SVaclav Hapla   if (!h->complexVal) {
1378c76551afSVaclav Hapla     H5T_class_t clazz = H5Tget_class(datatype);
1379c76551afSVaclav Hapla     if (clazz == H5T_FLOAT) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"File contains real numbers but PETSc is configured for complex. The conversion is not yet implemented. Configure with --with-scalar-type=real.");
13805a89fdf4SVaclav Hapla   }
13815a89fdf4SVaclav Hapla #else
13825a89fdf4SVaclav Hapla   if (h->complexVal) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"File contains complex numbers but PETSc not configured for them. Configure with --with-scalar-type=complex.");
13835a89fdf4SVaclav Hapla #endif
1384e9e90110SVaclav Hapla 
1385e9e90110SVaclav Hapla   ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr);
1386e9e90110SVaclav Hapla   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
1387e9e90110SVaclav Hapla   ierr = PetscViewerHDF5ReadSelectHyperslab_Private(viewer, h, map, &memspace);CHKERRQ(ierr);
1388e9e90110SVaclav Hapla 
13894fc17bcdSVaclav Hapla   unitsize = H5Tget_size(datatype);
13904fc17bcdSVaclav Hapla   if (h->complexVal) unitsize *= 2;
1391dff35581SVaclav Hapla   if (unitsize <= 0 || unitsize > PetscMax(sizeof(PetscInt),sizeof(PetscScalar))) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Sanity check failed: HDF5 function H5Tget_size(datatype) returned suspicious value %D",unitsize);
13924fc17bcdSVaclav Hapla   ierr = PetscMalloc(map->n*unitsize, &arr);CHKERRQ(ierr);
13934fc17bcdSVaclav Hapla 
1394eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadArray_Private(viewer, h, datatype, memspace, arr);CHKERRQ(ierr);
1395a25c73c6SVaclav Hapla   PetscStackCallHDF5(H5Sclose,(memspace));
1396eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr);
1397a25c73c6SVaclav Hapla   *newarr = arr;
1398a25c73c6SVaclav Hapla   PetscFunctionReturn(0);
1399a25c73c6SVaclav Hapla }
1400a25c73c6SVaclav Hapla 
1401c1aaad9cSVaclav Hapla /*@C
1402c1aaad9cSVaclav Hapla  PetscViewerHDF5ReadSizes - Read block size and global size of a vector (Vec or IS) stored in an HDF5 file.
1403c1aaad9cSVaclav Hapla 
1404c1aaad9cSVaclav Hapla   Input Parameters:
1405c1aaad9cSVaclav Hapla + viewer - The HDF5 viewer
1406c1aaad9cSVaclav Hapla - name   - The vector name
1407c1aaad9cSVaclav Hapla 
1408c1aaad9cSVaclav Hapla   Output Parameter:
1409c1aaad9cSVaclav Hapla + bs     - block size
1410c1aaad9cSVaclav Hapla - N      - global size
1411c1aaad9cSVaclav Hapla 
1412c1aaad9cSVaclav Hapla   Note:
1413c1aaad9cSVaclav Hapla   A vector is stored as an HDF5 dataspace with 1-4 dimensions in this order:
1414c1aaad9cSVaclav Hapla   1) # timesteps (optional), 2) # blocks, 3) # elements per block (optional), 4) real and imaginary part (only for complex).
1415c1aaad9cSVaclav Hapla 
1416c1aaad9cSVaclav Hapla   A vectors can be stored as a 2D dataspace even if its blocksize is 1; see PetscViewerHDF5SetBaseDimension2().
1417c1aaad9cSVaclav Hapla 
1418c1aaad9cSVaclav Hapla   Level: advanced
1419c1aaad9cSVaclav Hapla 
1420c1aaad9cSVaclav Hapla .seealso: PetscViewerHDF5Open(), VecLoad(), ISLoad(), VecGetSize(), ISGetSize(), PetscViewerHDF5SetBaseDimension2()
1421c1aaad9cSVaclav Hapla @*/
142269a06e7bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadSizes(PetscViewer viewer, const char name[], PetscInt *bs, PetscInt *N)
142369a06e7bSVaclav Hapla {
1424fbd37863SVaclav Hapla   HDF5ReadCtx    h=NULL;
14258374c777SVaclav Hapla   PetscLayout    map=NULL;
142669a06e7bSVaclav Hapla   PetscErrorCode ierr;
142769a06e7bSVaclav Hapla 
142869a06e7bSVaclav Hapla   PetscFunctionBegin;
1429c1aaad9cSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1430eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr);
1431eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr);
1432eb5a92b4SVaclav Hapla   ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr);
14338374c777SVaclav Hapla   if (bs) *bs = map->bs;
14348374c777SVaclav Hapla   if (N) *N = map->N;
14358374c777SVaclav Hapla   ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr);
1436a5e1feadSVaclav Hapla   PetscFunctionReturn(0);
1437a5e1feadSVaclav Hapla }
1438a5e1feadSVaclav Hapla 
1439a75e6a4aSMatthew G. Knepley /*
1440a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1441a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
1442a75e6a4aSMatthew G. Knepley */
1443d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1444a75e6a4aSMatthew G. Knepley 
1445a75e6a4aSMatthew G. Knepley /*@C
1446a75e6a4aSMatthew G. Knepley   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
1447a75e6a4aSMatthew G. Knepley 
1448a75e6a4aSMatthew G. Knepley   Collective on MPI_Comm
1449a75e6a4aSMatthew G. Knepley 
1450a75e6a4aSMatthew G. Knepley   Input Parameter:
1451a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer
1452a75e6a4aSMatthew G. Knepley 
1453a75e6a4aSMatthew G. Knepley   Level: intermediate
1454a75e6a4aSMatthew G. Knepley 
1455a75e6a4aSMatthew G. Knepley   Options Database Keys:
1456a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name>
1457a75e6a4aSMatthew G. Knepley 
1458a75e6a4aSMatthew G. Knepley   Environmental variables:
1459a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME
1460a75e6a4aSMatthew G. Knepley 
1461a75e6a4aSMatthew G. Knepley   Notes:
1462a75e6a4aSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
1463a75e6a4aSMatthew G. Knepley   an error code.  The HDF5 PetscViewer is usually used in the form
1464a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1465a75e6a4aSMatthew G. Knepley 
1466a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1467a75e6a4aSMatthew G. Knepley @*/
1468a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1469a75e6a4aSMatthew G. Knepley {
1470a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
1471a75e6a4aSMatthew G. Knepley   PetscBool      flg;
1472a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
1473a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
1474a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
1475a75e6a4aSMatthew G. Knepley 
1476a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
1477a75e6a4aSMatthew 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);}
1478a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
147912801b39SBarry Smith     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
1480a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1481a75e6a4aSMatthew G. Knepley   }
148247435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1483a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1484a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
1485a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
1486a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1487a75e6a4aSMatthew G. Knepley     if (!flg) {
1488a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname,"output.h5");
1489a75e6a4aSMatthew G. Knepley       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1490a75e6a4aSMatthew G. Knepley     }
1491a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
1492a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1493a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
1494a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
149547435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
1496a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1497a75e6a4aSMatthew G. Knepley   }
1498a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
1499a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1500a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
1501a75e6a4aSMatthew G. Knepley }
1502