xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 9a0502c650a1f9e8ff1ace346fad63beb9fe6c5d)
1af0996ceSBarry Smith #include <petsc/private/viewerimpl.h>    /*I   "petscsys.h"   I*/
2d70abbfaSBarry Smith #include <petscviewerhdf5.h>    /*I   "petscviewerhdf5.h"   I*/
35c6c1daeSBarry Smith 
45c6c1daeSBarry Smith typedef struct GroupList {
55c6c1daeSBarry Smith   const char       *name;
65c6c1daeSBarry Smith   struct GroupList *next;
75c6c1daeSBarry Smith } GroupList;
85c6c1daeSBarry Smith 
95c6c1daeSBarry Smith typedef struct {
105c6c1daeSBarry Smith   char          *filename;
115c6c1daeSBarry Smith   PetscFileMode btype;
125c6c1daeSBarry Smith   hid_t         file_id;
135c6c1daeSBarry Smith   PetscInt      timestep;
145c6c1daeSBarry Smith   GroupList     *groups;
1582ea9b62SBarry Smith   PetscBool     basedimension2;  /* save vectors and DMDA vectors with a dimension of at least 2 even if the bs/dof is 1 */
16*9a0502c6SHåkon Strandenes   PetscBool     spoutput;  /* write data in single precision even if PETSc is compiled with double precision PetscReal */
175c6c1daeSBarry Smith } PetscViewer_HDF5;
185c6c1daeSBarry Smith 
195c6c1daeSBarry Smith #undef __FUNCT__
2082ea9b62SBarry Smith #define __FUNCT__ "PetscViewerSetFromOptions_HDF5"
2182ea9b62SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptions *PetscOptionsObject,PetscViewer v)
2282ea9b62SBarry Smith {
2382ea9b62SBarry Smith   PetscErrorCode   ierr;
2482ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;
2582ea9b62SBarry Smith 
2682ea9b62SBarry Smith   PetscFunctionBegin;
2782ea9b62SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr);
2882ea9b62SBarry Smith   ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr);
29*9a0502c6SHåkon Strandenes   ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr);
3082ea9b62SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3182ea9b62SBarry Smith   PetscFunctionReturn(0);
3282ea9b62SBarry Smith }
3382ea9b62SBarry Smith 
3482ea9b62SBarry Smith #undef __FUNCT__
355c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerFileClose_HDF5"
365c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
375c6c1daeSBarry Smith {
385c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
395c6c1daeSBarry Smith   PetscErrorCode   ierr;
405c6c1daeSBarry Smith 
415c6c1daeSBarry Smith   PetscFunctionBegin;
425c6c1daeSBarry Smith   ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);
43729ab6d9SBarry Smith   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
445c6c1daeSBarry Smith   PetscFunctionReturn(0);
455c6c1daeSBarry Smith }
465c6c1daeSBarry Smith 
475c6c1daeSBarry Smith #undef __FUNCT__
485c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerDestroy_HDF5"
495c6c1daeSBarry Smith PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
505c6c1daeSBarry Smith {
515c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
525c6c1daeSBarry Smith   PetscErrorCode   ierr;
535c6c1daeSBarry Smith 
545c6c1daeSBarry Smith   PetscFunctionBegin;
555c6c1daeSBarry Smith   ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr);
565c6c1daeSBarry Smith   while (hdf5->groups) {
575c6c1daeSBarry Smith     GroupList *tmp = hdf5->groups->next;
585c6c1daeSBarry Smith 
595c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups->name);CHKERRQ(ierr);
605c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups);CHKERRQ(ierr);
615c6c1daeSBarry Smith     hdf5->groups = tmp;
625c6c1daeSBarry Smith   }
635c6c1daeSBarry Smith   ierr = PetscFree(hdf5);CHKERRQ(ierr);
640b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
650b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
665c6c1daeSBarry Smith   PetscFunctionReturn(0);
675c6c1daeSBarry Smith }
685c6c1daeSBarry Smith 
695c6c1daeSBarry Smith #undef __FUNCT__
705c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerFileSetMode_HDF5"
715c6c1daeSBarry Smith PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
725c6c1daeSBarry Smith {
735c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
745c6c1daeSBarry Smith 
755c6c1daeSBarry Smith   PetscFunctionBegin;
765c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
775c6c1daeSBarry Smith   hdf5->btype = type;
785c6c1daeSBarry Smith   PetscFunctionReturn(0);
795c6c1daeSBarry Smith }
805c6c1daeSBarry Smith 
815c6c1daeSBarry Smith #undef __FUNCT__
8282ea9b62SBarry Smith #define __FUNCT__ "PetscViewerHDF5SetBaseDimension2_HDF5"
8382ea9b62SBarry Smith PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
8482ea9b62SBarry Smith {
8582ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
8682ea9b62SBarry Smith 
8782ea9b62SBarry Smith   PetscFunctionBegin;
8882ea9b62SBarry Smith   hdf5->basedimension2 = flg;
8982ea9b62SBarry Smith   PetscFunctionReturn(0);
9082ea9b62SBarry Smith }
9182ea9b62SBarry Smith 
9282ea9b62SBarry Smith #undef __FUNCT__
9382ea9b62SBarry Smith #define __FUNCT__ "PetscViewerHDF5SetBaseDimension2"
9482ea9b62SBarry Smith /*@C
9582ea9b62SBarry Smith      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
9682ea9b62SBarry Smith        dimension of 2.
9782ea9b62SBarry Smith 
9882ea9b62SBarry Smith     Logically Collective on PetscViewer
9982ea9b62SBarry Smith 
10082ea9b62SBarry Smith   Input Parameters:
10182ea9b62SBarry Smith +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
10282ea9b62SBarry 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
10382ea9b62SBarry Smith 
10482ea9b62SBarry Smith   Options Database:
10582ea9b62SBarry 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
10682ea9b62SBarry Smith 
10782ea9b62SBarry Smith 
10882ea9b62SBarry Smith   Notes: Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
10982ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
11082ea9b62SBarry Smith 
11182ea9b62SBarry Smith   Level: intermediate
11282ea9b62SBarry Smith 
11382ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
11482ea9b62SBarry Smith 
11582ea9b62SBarry Smith @*/
11682ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
11782ea9b62SBarry Smith {
11882ea9b62SBarry Smith   PetscErrorCode ierr;
11982ea9b62SBarry Smith 
12082ea9b62SBarry Smith   PetscFunctionBegin;
12182ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
12282ea9b62SBarry Smith   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
12382ea9b62SBarry Smith   PetscFunctionReturn(0);
12482ea9b62SBarry Smith }
12582ea9b62SBarry Smith 
12682ea9b62SBarry Smith #undef __FUNCT__
12782ea9b62SBarry Smith #define __FUNCT__ "PetscViewerHDF5GetBaseDimension2"
12882ea9b62SBarry Smith /*@C
12982ea9b62SBarry Smith      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
13082ea9b62SBarry Smith        dimension of 2.
13182ea9b62SBarry Smith 
13282ea9b62SBarry Smith     Logically Collective on PetscViewer
13382ea9b62SBarry Smith 
13482ea9b62SBarry Smith   Input Parameter:
13582ea9b62SBarry Smith .  viewer - the PetscViewer, must be of type HDF5
13682ea9b62SBarry Smith 
13782ea9b62SBarry Smith   Output Parameter:
13882ea9b62SBarry 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
13982ea9b62SBarry Smith 
14082ea9b62SBarry Smith   Notes: Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
14182ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
14282ea9b62SBarry Smith 
14382ea9b62SBarry Smith   Level: intermediate
14482ea9b62SBarry Smith 
14582ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
14682ea9b62SBarry Smith 
14782ea9b62SBarry Smith @*/
14882ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
14982ea9b62SBarry Smith {
15082ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
15182ea9b62SBarry Smith 
15282ea9b62SBarry Smith   PetscFunctionBegin;
15382ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
15482ea9b62SBarry Smith   *flg = hdf5->basedimension2;
15582ea9b62SBarry Smith   PetscFunctionReturn(0);
15682ea9b62SBarry Smith }
15782ea9b62SBarry Smith 
15882ea9b62SBarry Smith #undef __FUNCT__
159*9a0502c6SHåkon Strandenes #define __FUNCT__ "PetscViewerHDF5SetSPOutput_HDF5"
160*9a0502c6SHåkon Strandenes PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
161*9a0502c6SHåkon Strandenes {
162*9a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
163*9a0502c6SHåkon Strandenes 
164*9a0502c6SHåkon Strandenes   PetscFunctionBegin;
165*9a0502c6SHåkon Strandenes   hdf5->spoutput = flg;
166*9a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
167*9a0502c6SHåkon Strandenes }
168*9a0502c6SHåkon Strandenes 
169*9a0502c6SHåkon Strandenes #undef __FUNCT__
170*9a0502c6SHåkon Strandenes #define __FUNCT__ "PetscViewerHDF5SetSPOutput"
171*9a0502c6SHåkon Strandenes /*@C
172*9a0502c6SHåkon Strandenes      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
173*9a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
174*9a0502c6SHåkon Strandenes 
175*9a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
176*9a0502c6SHåkon Strandenes 
177*9a0502c6SHåkon Strandenes   Input Parameters:
178*9a0502c6SHåkon Strandenes +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
179*9a0502c6SHåkon Strandenes -  flg - if PETSC_TRUE the data will be written to disk with single precision
180*9a0502c6SHåkon Strandenes 
181*9a0502c6SHåkon Strandenes   Options Database:
182*9a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
183*9a0502c6SHåkon Strandenes 
184*9a0502c6SHåkon Strandenes 
185*9a0502c6SHåkon Strandenes   Notes: Setting this option does not make any difference if PETSc is compiled with single precision
186*9a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
187*9a0502c6SHåkon Strandenes 
188*9a0502c6SHåkon Strandenes   Level: intermediate
189*9a0502c6SHåkon Strandenes 
190*9a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
191*9a0502c6SHåkon Strandenes           PetscReal
192*9a0502c6SHåkon Strandenes 
193*9a0502c6SHåkon Strandenes @*/
194*9a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
195*9a0502c6SHåkon Strandenes {
196*9a0502c6SHåkon Strandenes   PetscErrorCode ierr;
197*9a0502c6SHåkon Strandenes 
198*9a0502c6SHåkon Strandenes   PetscFunctionBegin;
199*9a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
200*9a0502c6SHåkon Strandenes   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
201*9a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
202*9a0502c6SHåkon Strandenes }
203*9a0502c6SHåkon Strandenes 
204*9a0502c6SHåkon Strandenes #undef __FUNCT__
205*9a0502c6SHåkon Strandenes #define __FUNCT__ "PetscViewerHDF5GetSPOutput"
206*9a0502c6SHåkon Strandenes /*@C
207*9a0502c6SHåkon Strandenes      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
208*9a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
209*9a0502c6SHåkon Strandenes 
210*9a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
211*9a0502c6SHåkon Strandenes 
212*9a0502c6SHåkon Strandenes   Input Parameter:
213*9a0502c6SHåkon Strandenes .  viewer - the PetscViewer, must be of type HDF5
214*9a0502c6SHåkon Strandenes 
215*9a0502c6SHåkon Strandenes   Output Parameter:
216*9a0502c6SHåkon Strandenes .  flg - if PETSC_TRUE the data will be written to disk with single precision
217*9a0502c6SHåkon Strandenes 
218*9a0502c6SHåkon Strandenes   Notes: Setting this option does not make any difference if PETSc is compiled with single precision
219*9a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
220*9a0502c6SHåkon Strandenes 
221*9a0502c6SHåkon Strandenes   Level: intermediate
222*9a0502c6SHåkon Strandenes 
223*9a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
224*9a0502c6SHåkon Strandenes           PetscReal
225*9a0502c6SHåkon Strandenes 
226*9a0502c6SHåkon Strandenes @*/
227*9a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
228*9a0502c6SHåkon Strandenes {
229*9a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
230*9a0502c6SHåkon Strandenes 
231*9a0502c6SHåkon Strandenes   PetscFunctionBegin;
232*9a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
233*9a0502c6SHåkon Strandenes   *flg = hdf5->spoutput;
234*9a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
235*9a0502c6SHåkon Strandenes }
236*9a0502c6SHåkon Strandenes 
237*9a0502c6SHåkon Strandenes #undef __FUNCT__
2385c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerFileSetName_HDF5"
2395c6c1daeSBarry Smith PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
2405c6c1daeSBarry Smith {
2415c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2425c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
2435c6c1daeSBarry Smith   MPI_Info          info = MPI_INFO_NULL;
2445c6c1daeSBarry Smith #endif
2455c6c1daeSBarry Smith   hid_t             plist_id;
2465c6c1daeSBarry Smith   PetscErrorCode    ierr;
2475c6c1daeSBarry Smith 
2485c6c1daeSBarry Smith   PetscFunctionBegin;
2495c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
2505c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
251729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
2525c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
253729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info));
2545c6c1daeSBarry Smith #endif
2555c6c1daeSBarry Smith   /* Create or open the file collectively */
2565c6c1daeSBarry Smith   switch (hdf5->btype) {
2575c6c1daeSBarry Smith   case FILE_MODE_READ:
258729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
2595c6c1daeSBarry Smith     break;
2605c6c1daeSBarry Smith   case FILE_MODE_APPEND:
261729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
2625c6c1daeSBarry Smith     break;
2635c6c1daeSBarry Smith   case FILE_MODE_WRITE:
264729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
2655c6c1daeSBarry Smith     break;
2665c6c1daeSBarry Smith   default:
2675c6c1daeSBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
2685c6c1daeSBarry Smith   }
2695c6c1daeSBarry Smith   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
270729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
2715c6c1daeSBarry Smith   PetscFunctionReturn(0);
2725c6c1daeSBarry Smith }
2735c6c1daeSBarry Smith 
2745c6c1daeSBarry Smith #undef __FUNCT__
2755c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerCreate_HDF5"
2768cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
2775c6c1daeSBarry Smith {
2785c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
2795c6c1daeSBarry Smith   PetscErrorCode   ierr;
2805c6c1daeSBarry Smith 
2815c6c1daeSBarry Smith   PetscFunctionBegin;
282b00a9115SJed Brown   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
2835c6c1daeSBarry Smith 
2845c6c1daeSBarry Smith   v->data                = (void*) hdf5;
2855c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
28682ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
2875c6c1daeSBarry Smith   v->ops->flush          = 0;
2885c6c1daeSBarry Smith   hdf5->btype            = (PetscFileMode) -1;
2895c6c1daeSBarry Smith   hdf5->filename         = 0;
2905c6c1daeSBarry Smith   hdf5->timestep         = -1;
2910298fd71SBarry Smith   hdf5->groups           = NULL;
2925c6c1daeSBarry Smith 
2930b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
2940b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
29582ea9b62SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
296*9a0502c6SHåkon Strandenes   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr);
2975c6c1daeSBarry Smith   PetscFunctionReturn(0);
2985c6c1daeSBarry Smith }
2995c6c1daeSBarry Smith 
3005c6c1daeSBarry Smith #undef __FUNCT__
3015c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5Open"
3025c6c1daeSBarry Smith /*@C
3035c6c1daeSBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
3045c6c1daeSBarry Smith 
3055c6c1daeSBarry Smith    Collective on MPI_Comm
3065c6c1daeSBarry Smith 
3075c6c1daeSBarry Smith    Input Parameters:
3085c6c1daeSBarry Smith +  comm - MPI communicator
3095c6c1daeSBarry Smith .  name - name of file
3105c6c1daeSBarry Smith -  type - type of file
3115c6c1daeSBarry Smith $    FILE_MODE_WRITE - create new file for binary output
3125c6c1daeSBarry Smith $    FILE_MODE_READ - open existing file for binary input
3135c6c1daeSBarry Smith $    FILE_MODE_APPEND - open existing file for binary output
3145c6c1daeSBarry Smith 
3155c6c1daeSBarry Smith    Output Parameter:
3165c6c1daeSBarry Smith .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
3175c6c1daeSBarry Smith 
31882ea9b62SBarry Smith   Options Database:
31982ea9b62SBarry 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
320*9a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
32182ea9b62SBarry Smith 
3225c6c1daeSBarry Smith    Level: beginner
3235c6c1daeSBarry Smith 
3245c6c1daeSBarry Smith    Note:
3255c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
3265c6c1daeSBarry Smith 
3275c6c1daeSBarry Smith    Concepts: HDF5 files
3285c6c1daeSBarry Smith    Concepts: PetscViewerHDF5^creating
3295c6c1daeSBarry Smith 
330*9a0502c6SHåkon Strandenes .seealso: PetscViewerASCIIOpen(), PetscViewerSetFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
331*9a0502c6SHåkon Strandenes           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
332*9a0502c6SHåkon Strandenes           MatLoad(), PetscFileMode, PetscViewer
3335c6c1daeSBarry Smith @*/
3345c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
3355c6c1daeSBarry Smith {
3365c6c1daeSBarry Smith   PetscErrorCode ierr;
3375c6c1daeSBarry Smith 
3385c6c1daeSBarry Smith   PetscFunctionBegin;
3395c6c1daeSBarry Smith   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
3405c6c1daeSBarry Smith   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
3415c6c1daeSBarry Smith   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
3425c6c1daeSBarry Smith   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
3435c6c1daeSBarry Smith   PetscFunctionReturn(0);
3445c6c1daeSBarry Smith }
3455c6c1daeSBarry Smith 
3465c6c1daeSBarry Smith #undef __FUNCT__
3475c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5GetFileId"
3485c6c1daeSBarry Smith /*@C
3495c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
3505c6c1daeSBarry Smith 
3515c6c1daeSBarry Smith   Not collective
3525c6c1daeSBarry Smith 
3535c6c1daeSBarry Smith   Input Parameter:
3545c6c1daeSBarry Smith . viewer - the PetscViewer
3555c6c1daeSBarry Smith 
3565c6c1daeSBarry Smith   Output Parameter:
3575c6c1daeSBarry Smith . file_id - The file id
3585c6c1daeSBarry Smith 
3595c6c1daeSBarry Smith   Level: intermediate
3605c6c1daeSBarry Smith 
3615c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open()
3625c6c1daeSBarry Smith @*/
3635c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
3645c6c1daeSBarry Smith {
3655c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3665c6c1daeSBarry Smith 
3675c6c1daeSBarry Smith   PetscFunctionBegin;
3685c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
3695c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
3705c6c1daeSBarry Smith   PetscFunctionReturn(0);
3715c6c1daeSBarry Smith }
3725c6c1daeSBarry Smith 
3735c6c1daeSBarry Smith #undef __FUNCT__
3745c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5PushGroup"
3755c6c1daeSBarry Smith /*@C
3765c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
3775c6c1daeSBarry Smith 
3785c6c1daeSBarry Smith   Not collective
3795c6c1daeSBarry Smith 
3805c6c1daeSBarry Smith   Input Parameters:
3815c6c1daeSBarry Smith + viewer - the PetscViewer
3825c6c1daeSBarry Smith - name - The group name
3835c6c1daeSBarry Smith 
3845c6c1daeSBarry Smith   Level: intermediate
3855c6c1daeSBarry Smith 
3865c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
3875c6c1daeSBarry Smith @*/
3885c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
3895c6c1daeSBarry Smith {
3905c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3915c6c1daeSBarry Smith   GroupList        *groupNode;
3925c6c1daeSBarry Smith   PetscErrorCode   ierr;
3935c6c1daeSBarry Smith 
3945c6c1daeSBarry Smith   PetscFunctionBegin;
3955c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
3965c6c1daeSBarry Smith   PetscValidCharPointer(name,2);
3975c6c1daeSBarry Smith   ierr = PetscMalloc(sizeof(GroupList), &groupNode);CHKERRQ(ierr);
3985c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
399a297a907SKarl Rupp 
4005c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
4015c6c1daeSBarry Smith   hdf5->groups    = groupNode;
4025c6c1daeSBarry Smith   PetscFunctionReturn(0);
4035c6c1daeSBarry Smith }
4045c6c1daeSBarry Smith 
4055c6c1daeSBarry Smith #undef __FUNCT__
4065c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5PopGroup"
4073ef9c667SSatish Balay /*@
4085c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
4095c6c1daeSBarry Smith 
4105c6c1daeSBarry Smith   Not collective
4115c6c1daeSBarry Smith 
4125c6c1daeSBarry Smith   Input Parameter:
4135c6c1daeSBarry Smith . viewer - the PetscViewer
4145c6c1daeSBarry Smith 
4155c6c1daeSBarry Smith   Level: intermediate
4165c6c1daeSBarry Smith 
4175c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup()
4185c6c1daeSBarry Smith @*/
4195c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
4205c6c1daeSBarry Smith {
4215c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4225c6c1daeSBarry Smith   GroupList        *groupNode;
4235c6c1daeSBarry Smith   PetscErrorCode   ierr;
4245c6c1daeSBarry Smith 
4255c6c1daeSBarry Smith   PetscFunctionBegin;
4265c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
42782f516ccSBarry Smith   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
4285c6c1daeSBarry Smith   groupNode    = hdf5->groups;
4295c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
4305c6c1daeSBarry Smith   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
4315c6c1daeSBarry Smith   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
4325c6c1daeSBarry Smith   PetscFunctionReturn(0);
4335c6c1daeSBarry Smith }
4345c6c1daeSBarry Smith 
4355c6c1daeSBarry Smith #undef __FUNCT__
4365c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5GetGroup"
4375c6c1daeSBarry Smith /*@C
4380298fd71SBarry Smith   PetscViewerHDF5GetGroup - Get the current HDF5 group for output. If none has been assigned, returns NULL.
4395c6c1daeSBarry Smith 
4405c6c1daeSBarry Smith   Not collective
4415c6c1daeSBarry Smith 
4425c6c1daeSBarry Smith   Input Parameter:
4435c6c1daeSBarry Smith . viewer - the PetscViewer
4445c6c1daeSBarry Smith 
4455c6c1daeSBarry Smith   Output Parameter:
4465c6c1daeSBarry Smith . name - The group name
4475c6c1daeSBarry Smith 
4485c6c1daeSBarry Smith   Level: intermediate
4495c6c1daeSBarry Smith 
4505c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup()
4515c6c1daeSBarry Smith @*/
4525c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
4535c6c1daeSBarry Smith {
4545c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
4555c6c1daeSBarry Smith 
4565c6c1daeSBarry Smith   PetscFunctionBegin;
4575c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
458c959eef4SJed Brown   PetscValidPointer(name,2);
459a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
4600298fd71SBarry Smith   else *name = NULL;
4615c6c1daeSBarry Smith   PetscFunctionReturn(0);
4625c6c1daeSBarry Smith }
4635c6c1daeSBarry Smith 
4645c6c1daeSBarry Smith #undef __FUNCT__
4655c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5IncrementTimestep"
4663ef9c667SSatish Balay /*@
4675c6c1daeSBarry Smith   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
4685c6c1daeSBarry Smith 
4695c6c1daeSBarry Smith   Not collective
4705c6c1daeSBarry Smith 
4715c6c1daeSBarry Smith   Input Parameter:
4725c6c1daeSBarry Smith . viewer - the PetscViewer
4735c6c1daeSBarry Smith 
4745c6c1daeSBarry Smith   Level: intermediate
4755c6c1daeSBarry Smith 
4765c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
4775c6c1daeSBarry Smith @*/
4785c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
4795c6c1daeSBarry Smith {
4805c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4815c6c1daeSBarry Smith 
4825c6c1daeSBarry Smith   PetscFunctionBegin;
4835c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
4845c6c1daeSBarry Smith   ++hdf5->timestep;
4855c6c1daeSBarry Smith   PetscFunctionReturn(0);
4865c6c1daeSBarry Smith }
4875c6c1daeSBarry Smith 
4885c6c1daeSBarry Smith #undef __FUNCT__
4895c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5SetTimestep"
4903ef9c667SSatish Balay /*@
4915c6c1daeSBarry Smith   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
4925c6c1daeSBarry Smith   of -1 disables blocking with timesteps.
4935c6c1daeSBarry Smith 
4945c6c1daeSBarry Smith   Not collective
4955c6c1daeSBarry Smith 
4965c6c1daeSBarry Smith   Input Parameters:
4975c6c1daeSBarry Smith + viewer - the PetscViewer
4985c6c1daeSBarry Smith - timestep - The timestep number
4995c6c1daeSBarry Smith 
5005c6c1daeSBarry Smith   Level: intermediate
5015c6c1daeSBarry Smith 
5025c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
5035c6c1daeSBarry Smith @*/
5045c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
5055c6c1daeSBarry Smith {
5065c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5075c6c1daeSBarry Smith 
5085c6c1daeSBarry Smith   PetscFunctionBegin;
5095c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5105c6c1daeSBarry Smith   hdf5->timestep = timestep;
5115c6c1daeSBarry Smith   PetscFunctionReturn(0);
5125c6c1daeSBarry Smith }
5135c6c1daeSBarry Smith 
5145c6c1daeSBarry Smith #undef __FUNCT__
5155c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5GetTimestep"
5163ef9c667SSatish Balay /*@
5175c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
5185c6c1daeSBarry Smith 
5195c6c1daeSBarry Smith   Not collective
5205c6c1daeSBarry Smith 
5215c6c1daeSBarry Smith   Input Parameter:
5225c6c1daeSBarry Smith . viewer - the PetscViewer
5235c6c1daeSBarry Smith 
5245c6c1daeSBarry Smith   Output Parameter:
5255c6c1daeSBarry Smith . timestep - The timestep number
5265c6c1daeSBarry Smith 
5275c6c1daeSBarry Smith   Level: intermediate
5285c6c1daeSBarry Smith 
5295c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
5305c6c1daeSBarry Smith @*/
5315c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
5325c6c1daeSBarry Smith {
5335c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5345c6c1daeSBarry Smith 
5355c6c1daeSBarry Smith   PetscFunctionBegin;
5365c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5375c6c1daeSBarry Smith   PetscValidPointer(timestep,2);
5385c6c1daeSBarry Smith   *timestep = hdf5->timestep;
5395c6c1daeSBarry Smith   PetscFunctionReturn(0);
5405c6c1daeSBarry Smith }
5415c6c1daeSBarry Smith 
54236bce6e8SMatthew G. Knepley #undef __FUNCT__
54336bce6e8SMatthew G. Knepley #define __FUNCT__ "PetscDataTypeToHDF5DataType"
54436bce6e8SMatthew G. Knepley /*@C
54536bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
54636bce6e8SMatthew G. Knepley 
54736bce6e8SMatthew G. Knepley   Not collective
54836bce6e8SMatthew G. Knepley 
54936bce6e8SMatthew G. Knepley   Input Parameter:
55036bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
55136bce6e8SMatthew G. Knepley 
55236bce6e8SMatthew G. Knepley   Output Parameter:
55336bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
55436bce6e8SMatthew G. Knepley 
55536bce6e8SMatthew G. Knepley   Level: advanced
55636bce6e8SMatthew G. Knepley 
55736bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
55836bce6e8SMatthew G. Knepley @*/
55936bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
56036bce6e8SMatthew G. Knepley {
56136bce6e8SMatthew G. Knepley   PetscFunctionBegin;
56236bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
56336bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
56436bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_LLONG;
56536bce6e8SMatthew G. Knepley #else
56636bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_INT;
56736bce6e8SMatthew G. Knepley #endif
56836bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
56936bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
57036bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
57136bce6e8SMatthew G. Knepley   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_DOUBLE;
57236bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_DOUBLE;
57336bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
57436bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
57536bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
5767e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
57736bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
57836bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
57936bce6e8SMatthew G. Knepley }
58036bce6e8SMatthew G. Knepley 
58136bce6e8SMatthew G. Knepley #undef __FUNCT__
58236bce6e8SMatthew G. Knepley #define __FUNCT__ "PetscHDF5DataTypeToPetscDataType"
58336bce6e8SMatthew G. Knepley /*@C
58436bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
58536bce6e8SMatthew G. Knepley 
58636bce6e8SMatthew G. Knepley   Not collective
58736bce6e8SMatthew G. Knepley 
58836bce6e8SMatthew G. Knepley   Input Parameter:
58936bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
59036bce6e8SMatthew G. Knepley 
59136bce6e8SMatthew G. Knepley   Output Parameter:
59236bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
59336bce6e8SMatthew G. Knepley 
59436bce6e8SMatthew G. Knepley   Level: advanced
59536bce6e8SMatthew G. Knepley 
59636bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
59736bce6e8SMatthew G. Knepley @*/
59836bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
59936bce6e8SMatthew G. Knepley {
60036bce6e8SMatthew G. Knepley   PetscFunctionBegin;
60136bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
60236bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
60336bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
60436bce6e8SMatthew G. Knepley #else
60536bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
60636bce6e8SMatthew G. Knepley #endif
60736bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
60836bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
60936bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
61036bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
61136bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
61236bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
6137e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
61436bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
61536bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
61636bce6e8SMatthew G. Knepley }
61736bce6e8SMatthew G. Knepley 
61836bce6e8SMatthew G. Knepley #undef __FUNCT__
61936bce6e8SMatthew G. Knepley #define __FUNCT__ "PetscViewerHDF5WriteAttribute"
62036bce6e8SMatthew G. Knepley /*@
62136bce6e8SMatthew G. Knepley  PetscViewerHDF5WriteAttribute - Write a scalar attribute
62236bce6e8SMatthew G. Knepley 
62336bce6e8SMatthew G. Knepley   Input Parameters:
62436bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer
62536bce6e8SMatthew G. Knepley . parent - The parent name
62636bce6e8SMatthew G. Knepley . name   - The attribute name
62736bce6e8SMatthew G. Knepley . datatype - The attribute type
62836bce6e8SMatthew G. Knepley - value    - The attribute value
62936bce6e8SMatthew G. Knepley 
63036bce6e8SMatthew G. Knepley   Level: advanced
63136bce6e8SMatthew G. Knepley 
632e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute()
63336bce6e8SMatthew G. Knepley @*/
63436bce6e8SMatthew G. Knepley PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
63536bce6e8SMatthew G. Knepley {
636729ab6d9SBarry Smith   hid_t          h5, dataspace, dataset, attribute, dtype;
63736bce6e8SMatthew G. Knepley   PetscErrorCode ierr;
63836bce6e8SMatthew G. Knepley 
63936bce6e8SMatthew G. Knepley   PetscFunctionBegin;
6405cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
64136bce6e8SMatthew G. Knepley   PetscValidPointer(parent, 2);
64236bce6e8SMatthew G. Knepley   PetscValidPointer(name, 3);
64336bce6e8SMatthew G. Knepley   PetscValidPointer(value, 4);
64436bce6e8SMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
6457e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
6467e97c476SMatthew G. Knepley     size_t len;
6477e97c476SMatthew G. Knepley     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
648729ab6d9SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
6497e97c476SMatthew G. Knepley   }
65036bce6e8SMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
651729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
65236bce6e8SMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
653729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataset,H5Dopen2,(h5, parent, H5P_DEFAULT));
654729ab6d9SBarry Smith   PetscStackCallHDF5Return(attribute,H5Acreate2,(dataset, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
65536bce6e8SMatthew G. Knepley #else
656729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataset,H5Dopen,(h5, parent));
657729ab6d9SBarry Smith   PetscStackCallHDF5Return(attribute,H5Acreate,(dataset, name, dtype, dataspace, H5P_DEFAULT));
65836bce6e8SMatthew G. Knepley #endif
659729ab6d9SBarry Smith   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
660729ab6d9SBarry Smith   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
661729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
662729ab6d9SBarry Smith   PetscStackCallHDF5(H5Dclose,(dataset));
663729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(dataspace));
66436bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
66536bce6e8SMatthew G. Knepley }
66636bce6e8SMatthew G. Knepley 
667ad0c61feSMatthew G. Knepley #undef __FUNCT__
668ad0c61feSMatthew G. Knepley #define __FUNCT__ "PetscViewerHDF5ReadAttribute"
669ad0c61feSMatthew G. Knepley /*@
670ad0c61feSMatthew G. Knepley  PetscViewerHDF5ReadAttribute - Read a scalar attribute
671ad0c61feSMatthew G. Knepley 
672ad0c61feSMatthew G. Knepley   Input Parameters:
673ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer
674ad0c61feSMatthew G. Knepley . parent - The parent name
675ad0c61feSMatthew G. Knepley . name   - The attribute name
676ad0c61feSMatthew G. Knepley - datatype - The attribute type
677ad0c61feSMatthew G. Knepley 
678ad0c61feSMatthew G. Knepley   Output Parameter:
679ad0c61feSMatthew G. Knepley . value    - The attribute value
680ad0c61feSMatthew G. Knepley 
681ad0c61feSMatthew G. Knepley   Level: advanced
682ad0c61feSMatthew G. Knepley 
683e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute()
684ad0c61feSMatthew G. Knepley @*/
685ad0c61feSMatthew G. Knepley PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, void *value)
686ad0c61feSMatthew G. Knepley {
687729ab6d9SBarry Smith   hid_t          h5, dataspace, dataset, attribute, dtype;
688ad0c61feSMatthew G. Knepley   PetscErrorCode ierr;
689ad0c61feSMatthew G. Knepley 
690ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
6915cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
692ad0c61feSMatthew G. Knepley   PetscValidPointer(parent, 2);
693ad0c61feSMatthew G. Knepley   PetscValidPointer(name, 3);
694ad0c61feSMatthew G. Knepley   PetscValidPointer(value, 4);
695ad0c61feSMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
696ad0c61feSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
697729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
698ad0c61feSMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
699729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataset,H5Dopen2,(h5, parent, H5P_DEFAULT));
700ad0c61feSMatthew G. Knepley #else
701729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataset,H5Dopen,(h5, parent));
702ad0c61feSMatthew G. Knepley #endif
703729ab6d9SBarry Smith   PetscStackCallHDF5Return(attribute,H5Aopen_name,(dataset, name));
70470efba33SBarry Smith   PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
705729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
706729ab6d9SBarry Smith   PetscStackCallHDF5(H5Dclose,(dataset));
707729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(dataspace));
708ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
709ad0c61feSMatthew G. Knepley }
710ad0c61feSMatthew G. Knepley 
711e7bdbf8eSMatthew G. Knepley #undef __FUNCT__
7125cdeb108SMatthew G. Knepley #define __FUNCT__ "PetscViewerHDF5HasObject"
7135cdeb108SMatthew G. Knepley static PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, const char name[], H5O_type_t otype, PetscBool *has)
7145cdeb108SMatthew G. Knepley {
7155cdeb108SMatthew G. Knepley   hid_t          h5;
7165cdeb108SMatthew G. Knepley   PetscErrorCode ierr;
7175cdeb108SMatthew G. Knepley 
7185cdeb108SMatthew G. Knepley   PetscFunctionBegin;
7195cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7205cdeb108SMatthew G. Knepley   PetscValidPointer(name, 2);
7215cdeb108SMatthew G. Knepley   PetscValidPointer(has, 3);
7225cdeb108SMatthew G. Knepley   *has = PETSC_FALSE;
723c0cd0301SJed Brown   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
7245cdeb108SMatthew G. Knepley   if (H5Lexists(h5, name, H5P_DEFAULT)) {
7255cdeb108SMatthew G. Knepley     H5O_info_t info;
7265cdeb108SMatthew G. Knepley     hid_t      obj;
7275cdeb108SMatthew G. Knepley 
728729ab6d9SBarry Smith     PetscStackCallHDF5Return(obj,H5Oopen,(h5, name, H5P_DEFAULT));
729729ab6d9SBarry Smith     PetscStackCallHDF5(H5Oget_info,(obj, &info));
7305cdeb108SMatthew G. Knepley     if (otype == info.type) *has = PETSC_TRUE;
731729ab6d9SBarry Smith     PetscStackCallHDF5(H5Oclose,(obj));
7325cdeb108SMatthew G. Knepley   }
7335cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
7345cdeb108SMatthew G. Knepley }
7355cdeb108SMatthew G. Knepley 
7365cdeb108SMatthew G. Knepley #undef __FUNCT__
737e7bdbf8eSMatthew G. Knepley #define __FUNCT__ "PetscViewerHDF5HasAttribute"
738e7bdbf8eSMatthew G. Knepley /*@
739e7bdbf8eSMatthew G. Knepley  PetscViewerHDF5HasAttribute - Check whether a scalar attribute exists
740e7bdbf8eSMatthew G. Knepley 
741e7bdbf8eSMatthew G. Knepley   Input Parameters:
742e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer
743e7bdbf8eSMatthew G. Knepley . parent - The parent name
744e7bdbf8eSMatthew G. Knepley - name   - The attribute name
745e7bdbf8eSMatthew G. Knepley 
746e7bdbf8eSMatthew G. Knepley   Output Parameter:
747e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
748e7bdbf8eSMatthew G. Knepley 
749e7bdbf8eSMatthew G. Knepley   Level: advanced
750e7bdbf8eSMatthew G. Knepley 
751e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute()
752e7bdbf8eSMatthew G. Knepley @*/
753e7bdbf8eSMatthew G. Knepley PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
754e7bdbf8eSMatthew G. Knepley {
755729ab6d9SBarry Smith   hid_t          h5, dataset;
7565cdeb108SMatthew G. Knepley   htri_t         hhas;
7575cdeb108SMatthew G. Knepley   PetscBool      exists;
758e7bdbf8eSMatthew G. Knepley   PetscErrorCode ierr;
759e7bdbf8eSMatthew G. Knepley 
760e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
7615cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
762e7bdbf8eSMatthew G. Knepley   PetscValidPointer(parent, 2);
763e7bdbf8eSMatthew G. Knepley   PetscValidPointer(name, 3);
764e7bdbf8eSMatthew G. Knepley   PetscValidPointer(has, 4);
765e7bdbf8eSMatthew G. Knepley   *has = PETSC_FALSE;
766e7bdbf8eSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
7675cdeb108SMatthew G. Knepley   ierr = PetscViewerHDF5HasObject(viewer, parent, H5O_TYPE_DATASET, &exists);CHKERRQ(ierr);
7685cdeb108SMatthew G. Knepley   if (exists) {
769e7bdbf8eSMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
770729ab6d9SBarry Smith     PetscStackCall("H5Dopen2",dataset = H5Dopen2(h5, parent, H5P_DEFAULT));
771e7bdbf8eSMatthew G. Knepley #else
772729ab6d9SBarry Smith     PetscStackCall("H5Dopen",dataset = H5Dopen(h5, parent));
773e7bdbf8eSMatthew G. Knepley #endif
774729ab6d9SBarry Smith     if (dataset < 0) PetscFunctionReturn(0);
775729ab6d9SBarry Smith     PetscStackCall("H5Aexists",hhas = H5Aexists(dataset, name));
776729ab6d9SBarry Smith     if (hhas < 0) {
777729ab6d9SBarry Smith       PetscStackCallHDF5(H5Dclose,(dataset));
778729ab6d9SBarry Smith       PetscFunctionReturn(0);
779729ab6d9SBarry Smith     }
780729ab6d9SBarry Smith     PetscStackCallHDF5(H5Dclose,(dataset));
7815cdeb108SMatthew G. Knepley     *has = hhas ? PETSC_TRUE : PETSC_FALSE;
7825cdeb108SMatthew G. Knepley   }
783e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
784e7bdbf8eSMatthew G. Knepley }
785e7bdbf8eSMatthew G. Knepley 
786a75e6a4aSMatthew G. Knepley /*
787a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
788a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
789a75e6a4aSMatthew G. Knepley */
790a75e6a4aSMatthew G. Knepley static int Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
791a75e6a4aSMatthew G. Knepley 
792a75e6a4aSMatthew G. Knepley #undef __FUNCT__
793a75e6a4aSMatthew G. Knepley #define __FUNCT__ "PETSC_VIEWER_HDF5_"
794a75e6a4aSMatthew G. Knepley /*@C
795a75e6a4aSMatthew G. Knepley   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
796a75e6a4aSMatthew G. Knepley 
797a75e6a4aSMatthew G. Knepley   Collective on MPI_Comm
798a75e6a4aSMatthew G. Knepley 
799a75e6a4aSMatthew G. Knepley   Input Parameter:
800a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer
801a75e6a4aSMatthew G. Knepley 
802a75e6a4aSMatthew G. Knepley   Level: intermediate
803a75e6a4aSMatthew G. Knepley 
804a75e6a4aSMatthew G. Knepley   Options Database Keys:
805a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name>
806a75e6a4aSMatthew G. Knepley 
807a75e6a4aSMatthew G. Knepley   Environmental variables:
808a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME
809a75e6a4aSMatthew G. Knepley 
810a75e6a4aSMatthew G. Knepley   Notes:
811a75e6a4aSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
812a75e6a4aSMatthew G. Knepley   an error code.  The HDF5 PetscViewer is usually used in the form
813a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
814a75e6a4aSMatthew G. Knepley 
815a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
816a75e6a4aSMatthew G. Knepley @*/
817a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
818a75e6a4aSMatthew G. Knepley {
819a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
820a75e6a4aSMatthew G. Knepley   PetscBool      flg;
821a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
822a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
823a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
824a75e6a4aSMatthew G. Knepley 
825a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
826a75e6a4aSMatthew 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);}
827a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
828a75e6a4aSMatthew G. Knepley     ierr = MPI_Keyval_create(MPI_NULL_COPY_FN,MPI_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
829a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
830a75e6a4aSMatthew G. Knepley   }
831a75e6a4aSMatthew G. Knepley   ierr = MPI_Attr_get(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
832a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
833a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
834a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
835a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
836a75e6a4aSMatthew G. Knepley     if (!flg) {
837a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname,"output.h5");
838a75e6a4aSMatthew G. Knepley       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
839a75e6a4aSMatthew G. Knepley     }
840a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
841a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
842a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
843a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
844a75e6a4aSMatthew G. Knepley     ierr = MPI_Attr_put(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
845a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
846a75e6a4aSMatthew G. Knepley   }
847a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
848a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
849a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
850a75e6a4aSMatthew G. Knepley }
851