xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 95dccacacae8a8fc0b691f9b4fba69a249b61188)
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 */
169a0502c6SHå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"
214416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *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);
299a0502c6SHå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);
65d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr);
660b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
675c6c1daeSBarry Smith   PetscFunctionReturn(0);
685c6c1daeSBarry Smith }
695c6c1daeSBarry Smith 
705c6c1daeSBarry Smith #undef __FUNCT__
715c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerFileSetMode_HDF5"
725c6c1daeSBarry Smith PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
735c6c1daeSBarry Smith {
745c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
755c6c1daeSBarry Smith 
765c6c1daeSBarry Smith   PetscFunctionBegin;
775c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
785c6c1daeSBarry Smith   hdf5->btype = type;
795c6c1daeSBarry Smith   PetscFunctionReturn(0);
805c6c1daeSBarry Smith }
815c6c1daeSBarry Smith 
825c6c1daeSBarry Smith #undef __FUNCT__
8382ea9b62SBarry Smith #define __FUNCT__ "PetscViewerHDF5SetBaseDimension2_HDF5"
8482ea9b62SBarry Smith PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
8582ea9b62SBarry Smith {
8682ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
8782ea9b62SBarry Smith 
8882ea9b62SBarry Smith   PetscFunctionBegin;
8982ea9b62SBarry Smith   hdf5->basedimension2 = flg;
9082ea9b62SBarry Smith   PetscFunctionReturn(0);
9182ea9b62SBarry Smith }
9282ea9b62SBarry Smith 
9382ea9b62SBarry Smith #undef __FUNCT__
9482ea9b62SBarry Smith #define __FUNCT__ "PetscViewerHDF5SetBaseDimension2"
95df863907SAlex Fikl /*@
9682ea9b62SBarry Smith      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
9782ea9b62SBarry Smith        dimension of 2.
9882ea9b62SBarry Smith 
9982ea9b62SBarry Smith     Logically Collective on PetscViewer
10082ea9b62SBarry Smith 
10182ea9b62SBarry Smith   Input Parameters:
10282ea9b62SBarry Smith +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
10382ea9b62SBarry 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
10482ea9b62SBarry Smith 
10582ea9b62SBarry Smith   Options Database:
10682ea9b62SBarry 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
10782ea9b62SBarry Smith 
10882ea9b62SBarry Smith 
10982ea9b62SBarry 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
11082ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
11182ea9b62SBarry Smith 
11282ea9b62SBarry Smith   Level: intermediate
11382ea9b62SBarry Smith 
11482ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
11582ea9b62SBarry Smith 
11682ea9b62SBarry Smith @*/
11782ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
11882ea9b62SBarry Smith {
11982ea9b62SBarry Smith   PetscErrorCode ierr;
12082ea9b62SBarry Smith 
12182ea9b62SBarry Smith   PetscFunctionBegin;
12282ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
12382ea9b62SBarry Smith   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
12482ea9b62SBarry Smith   PetscFunctionReturn(0);
12582ea9b62SBarry Smith }
12682ea9b62SBarry Smith 
12782ea9b62SBarry Smith #undef __FUNCT__
12882ea9b62SBarry Smith #define __FUNCT__ "PetscViewerHDF5GetBaseDimension2"
129df863907SAlex Fikl /*@
13082ea9b62SBarry Smith      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
13182ea9b62SBarry Smith        dimension of 2.
13282ea9b62SBarry Smith 
13382ea9b62SBarry Smith     Logically Collective on PetscViewer
13482ea9b62SBarry Smith 
13582ea9b62SBarry Smith   Input Parameter:
13682ea9b62SBarry Smith .  viewer - the PetscViewer, must be of type HDF5
13782ea9b62SBarry Smith 
13882ea9b62SBarry Smith   Output Parameter:
13982ea9b62SBarry 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
14082ea9b62SBarry Smith 
14182ea9b62SBarry 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
14282ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
14382ea9b62SBarry Smith 
14482ea9b62SBarry Smith   Level: intermediate
14582ea9b62SBarry Smith 
14682ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
14782ea9b62SBarry Smith 
14882ea9b62SBarry Smith @*/
14982ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
15082ea9b62SBarry Smith {
15182ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
15282ea9b62SBarry Smith 
15382ea9b62SBarry Smith   PetscFunctionBegin;
15482ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
15582ea9b62SBarry Smith   *flg = hdf5->basedimension2;
15682ea9b62SBarry Smith   PetscFunctionReturn(0);
15782ea9b62SBarry Smith }
15882ea9b62SBarry Smith 
15982ea9b62SBarry Smith #undef __FUNCT__
1609a0502c6SHåkon Strandenes #define __FUNCT__ "PetscViewerHDF5SetSPOutput_HDF5"
1619a0502c6SHåkon Strandenes PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
1629a0502c6SHåkon Strandenes {
1639a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1649a0502c6SHåkon Strandenes 
1659a0502c6SHåkon Strandenes   PetscFunctionBegin;
1669a0502c6SHåkon Strandenes   hdf5->spoutput = flg;
1679a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
1689a0502c6SHåkon Strandenes }
1699a0502c6SHåkon Strandenes 
1709a0502c6SHåkon Strandenes #undef __FUNCT__
1719a0502c6SHåkon Strandenes #define __FUNCT__ "PetscViewerHDF5SetSPOutput"
172df863907SAlex Fikl /*@
1739a0502c6SHåkon Strandenes      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
1749a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
1759a0502c6SHåkon Strandenes 
1769a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
1779a0502c6SHåkon Strandenes 
1789a0502c6SHåkon Strandenes   Input Parameters:
1799a0502c6SHåkon Strandenes +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
1809a0502c6SHåkon Strandenes -  flg - if PETSC_TRUE the data will be written to disk with single precision
1819a0502c6SHåkon Strandenes 
1829a0502c6SHåkon Strandenes   Options Database:
1839a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
1849a0502c6SHåkon Strandenes 
1859a0502c6SHåkon Strandenes 
1869a0502c6SHåkon Strandenes   Notes: Setting this option does not make any difference if PETSc is compiled with single precision
1879a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
1889a0502c6SHåkon Strandenes 
1899a0502c6SHåkon Strandenes   Level: intermediate
1909a0502c6SHåkon Strandenes 
1919a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
1929a0502c6SHåkon Strandenes           PetscReal
1939a0502c6SHåkon Strandenes 
1949a0502c6SHåkon Strandenes @*/
1959a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
1969a0502c6SHåkon Strandenes {
1979a0502c6SHåkon Strandenes   PetscErrorCode ierr;
1989a0502c6SHåkon Strandenes 
1999a0502c6SHåkon Strandenes   PetscFunctionBegin;
2009a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2019a0502c6SHåkon Strandenes   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
2029a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2039a0502c6SHåkon Strandenes }
2049a0502c6SHåkon Strandenes 
2059a0502c6SHåkon Strandenes #undef __FUNCT__
2069a0502c6SHåkon Strandenes #define __FUNCT__ "PetscViewerHDF5GetSPOutput"
207df863907SAlex Fikl /*@
2089a0502c6SHåkon Strandenes      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
2099a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2109a0502c6SHåkon Strandenes 
2119a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2129a0502c6SHåkon Strandenes 
2139a0502c6SHåkon Strandenes   Input Parameter:
2149a0502c6SHåkon Strandenes .  viewer - the PetscViewer, must be of type HDF5
2159a0502c6SHåkon Strandenes 
2169a0502c6SHåkon Strandenes   Output Parameter:
2179a0502c6SHåkon Strandenes .  flg - if PETSC_TRUE the data will be written to disk with single precision
2189a0502c6SHåkon Strandenes 
2199a0502c6SHåkon Strandenes   Notes: Setting this option does not make any difference if PETSc is compiled with single precision
2209a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2219a0502c6SHåkon Strandenes 
2229a0502c6SHåkon Strandenes   Level: intermediate
2239a0502c6SHåkon Strandenes 
2249a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2259a0502c6SHåkon Strandenes           PetscReal
2269a0502c6SHåkon Strandenes 
2279a0502c6SHåkon Strandenes @*/
2289a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
2299a0502c6SHåkon Strandenes {
2309a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2319a0502c6SHåkon Strandenes 
2329a0502c6SHåkon Strandenes   PetscFunctionBegin;
2339a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2349a0502c6SHåkon Strandenes   *flg = hdf5->spoutput;
2359a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2369a0502c6SHåkon Strandenes }
2379a0502c6SHåkon Strandenes 
2389a0502c6SHåkon Strandenes #undef __FUNCT__
2395c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerFileSetName_HDF5"
2405c6c1daeSBarry Smith PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
2415c6c1daeSBarry Smith {
2425c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2435c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
2445c6c1daeSBarry Smith   MPI_Info          info = MPI_INFO_NULL;
2455c6c1daeSBarry Smith #endif
2465c6c1daeSBarry Smith   hid_t             plist_id;
2475c6c1daeSBarry Smith   PetscErrorCode    ierr;
2485c6c1daeSBarry Smith 
2495c6c1daeSBarry Smith   PetscFunctionBegin;
250f683cc58SToby Isaac   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
251f683cc58SToby Isaac   if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);}
2525c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
2535c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
254729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
2555c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
256729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info));
2575c6c1daeSBarry Smith #endif
2585c6c1daeSBarry Smith   /* Create or open the file collectively */
2595c6c1daeSBarry Smith   switch (hdf5->btype) {
2605c6c1daeSBarry Smith   case FILE_MODE_READ:
261729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
2625c6c1daeSBarry Smith     break;
2635c6c1daeSBarry Smith   case FILE_MODE_APPEND:
264729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
2655c6c1daeSBarry Smith     break;
2665c6c1daeSBarry Smith   case FILE_MODE_WRITE:
267729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
2685c6c1daeSBarry Smith     break;
2695c6c1daeSBarry Smith   default:
2705c6c1daeSBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
2715c6c1daeSBarry Smith   }
2725c6c1daeSBarry Smith   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
273729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
2745c6c1daeSBarry Smith   PetscFunctionReturn(0);
2755c6c1daeSBarry Smith }
2765c6c1daeSBarry Smith 
2775c6c1daeSBarry Smith #undef __FUNCT__
278d1232d7fSToby Isaac #define __FUNCT__ "PetscViewerFileGetName_HDF5"
279d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
280d1232d7fSToby Isaac {
281d1232d7fSToby Isaac   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
282d1232d7fSToby Isaac 
283d1232d7fSToby Isaac   PetscFunctionBegin;
284d1232d7fSToby Isaac   *name = vhdf5->filename;
285d1232d7fSToby Isaac   PetscFunctionReturn(0);
286d1232d7fSToby Isaac }
287d1232d7fSToby Isaac 
288d1232d7fSToby Isaac 
289d1232d7fSToby Isaac #undef __FUNCT__
2905c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerCreate_HDF5"
2918cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
2925c6c1daeSBarry Smith {
2935c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
2945c6c1daeSBarry Smith   PetscErrorCode   ierr;
2955c6c1daeSBarry Smith 
2965c6c1daeSBarry Smith   PetscFunctionBegin;
297b00a9115SJed Brown   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
2985c6c1daeSBarry Smith 
2995c6c1daeSBarry Smith   v->data                = (void*) hdf5;
3005c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
30182ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
3025c6c1daeSBarry Smith   v->ops->flush          = 0;
3035c6c1daeSBarry Smith   hdf5->btype            = (PetscFileMode) -1;
3045c6c1daeSBarry Smith   hdf5->filename         = 0;
3055c6c1daeSBarry Smith   hdf5->timestep         = -1;
3060298fd71SBarry Smith   hdf5->groups           = NULL;
3075c6c1daeSBarry Smith 
3080b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
309d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr);
3100b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
31182ea9b62SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
3129a0502c6SHåkon Strandenes   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr);
3135c6c1daeSBarry Smith   PetscFunctionReturn(0);
3145c6c1daeSBarry Smith }
3155c6c1daeSBarry Smith 
3165c6c1daeSBarry Smith #undef __FUNCT__
3175c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5Open"
3185c6c1daeSBarry Smith /*@C
3195c6c1daeSBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
3205c6c1daeSBarry Smith 
3215c6c1daeSBarry Smith    Collective on MPI_Comm
3225c6c1daeSBarry Smith 
3235c6c1daeSBarry Smith    Input Parameters:
3245c6c1daeSBarry Smith +  comm - MPI communicator
3255c6c1daeSBarry Smith .  name - name of file
3265c6c1daeSBarry Smith -  type - type of file
3275c6c1daeSBarry Smith $    FILE_MODE_WRITE - create new file for binary output
3285c6c1daeSBarry Smith $    FILE_MODE_READ - open existing file for binary input
3295c6c1daeSBarry Smith $    FILE_MODE_APPEND - open existing file for binary output
3305c6c1daeSBarry Smith 
3315c6c1daeSBarry Smith    Output Parameter:
3325c6c1daeSBarry Smith .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
3335c6c1daeSBarry Smith 
33482ea9b62SBarry Smith   Options Database:
33582ea9b62SBarry 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
3369a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
33782ea9b62SBarry Smith 
3385c6c1daeSBarry Smith    Level: beginner
3395c6c1daeSBarry Smith 
3405c6c1daeSBarry Smith    Note:
3415c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
3425c6c1daeSBarry Smith 
3435c6c1daeSBarry Smith    Concepts: HDF5 files
3445c6c1daeSBarry Smith    Concepts: PetscViewerHDF5^creating
3455c6c1daeSBarry Smith 
3466a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
3479a0502c6SHåkon Strandenes           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
3489a0502c6SHåkon Strandenes           MatLoad(), PetscFileMode, PetscViewer
3495c6c1daeSBarry Smith @*/
3505c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
3515c6c1daeSBarry Smith {
3525c6c1daeSBarry Smith   PetscErrorCode ierr;
3535c6c1daeSBarry Smith 
3545c6c1daeSBarry Smith   PetscFunctionBegin;
3555c6c1daeSBarry Smith   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
3565c6c1daeSBarry Smith   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
3575c6c1daeSBarry Smith   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
3585c6c1daeSBarry Smith   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
3595c6c1daeSBarry Smith   PetscFunctionReturn(0);
3605c6c1daeSBarry Smith }
3615c6c1daeSBarry Smith 
3625c6c1daeSBarry Smith #undef __FUNCT__
3635c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5GetFileId"
3645c6c1daeSBarry Smith /*@C
3655c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
3665c6c1daeSBarry Smith 
3675c6c1daeSBarry Smith   Not collective
3685c6c1daeSBarry Smith 
3695c6c1daeSBarry Smith   Input Parameter:
3705c6c1daeSBarry Smith . viewer - the PetscViewer
3715c6c1daeSBarry Smith 
3725c6c1daeSBarry Smith   Output Parameter:
3735c6c1daeSBarry Smith . file_id - The file id
3745c6c1daeSBarry Smith 
3755c6c1daeSBarry Smith   Level: intermediate
3765c6c1daeSBarry Smith 
3775c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open()
3785c6c1daeSBarry Smith @*/
3795c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
3805c6c1daeSBarry Smith {
3815c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3825c6c1daeSBarry Smith 
3835c6c1daeSBarry Smith   PetscFunctionBegin;
3845c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
3855c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
3865c6c1daeSBarry Smith   PetscFunctionReturn(0);
3875c6c1daeSBarry Smith }
3885c6c1daeSBarry Smith 
3895c6c1daeSBarry Smith #undef __FUNCT__
3905c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5PushGroup"
3915c6c1daeSBarry Smith /*@C
3925c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
3935c6c1daeSBarry Smith 
3945c6c1daeSBarry Smith   Not collective
3955c6c1daeSBarry Smith 
3965c6c1daeSBarry Smith   Input Parameters:
3975c6c1daeSBarry Smith + viewer - the PetscViewer
3985c6c1daeSBarry Smith - name - The group name
3995c6c1daeSBarry Smith 
4005c6c1daeSBarry Smith   Level: intermediate
4015c6c1daeSBarry Smith 
4025c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
4035c6c1daeSBarry Smith @*/
4045c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
4055c6c1daeSBarry Smith {
4065c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4075c6c1daeSBarry Smith   GroupList        *groupNode;
4085c6c1daeSBarry Smith   PetscErrorCode   ierr;
4095c6c1daeSBarry Smith 
4105c6c1daeSBarry Smith   PetscFunctionBegin;
4115c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
4125c6c1daeSBarry Smith   PetscValidCharPointer(name,2);
413*95dccacaSBarry Smith   ierr = PetscNew(&groupNode);CHKERRQ(ierr);
4145c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
415a297a907SKarl Rupp 
4165c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
4175c6c1daeSBarry Smith   hdf5->groups    = groupNode;
4185c6c1daeSBarry Smith   PetscFunctionReturn(0);
4195c6c1daeSBarry Smith }
4205c6c1daeSBarry Smith 
4215c6c1daeSBarry Smith #undef __FUNCT__
4225c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5PopGroup"
4233ef9c667SSatish Balay /*@
4245c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
4255c6c1daeSBarry Smith 
4265c6c1daeSBarry Smith   Not collective
4275c6c1daeSBarry Smith 
4285c6c1daeSBarry Smith   Input Parameter:
4295c6c1daeSBarry Smith . viewer - the PetscViewer
4305c6c1daeSBarry Smith 
4315c6c1daeSBarry Smith   Level: intermediate
4325c6c1daeSBarry Smith 
4335c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup()
4345c6c1daeSBarry Smith @*/
4355c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
4365c6c1daeSBarry Smith {
4375c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4385c6c1daeSBarry Smith   GroupList        *groupNode;
4395c6c1daeSBarry Smith   PetscErrorCode   ierr;
4405c6c1daeSBarry Smith 
4415c6c1daeSBarry Smith   PetscFunctionBegin;
4425c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
44382f516ccSBarry Smith   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
4445c6c1daeSBarry Smith   groupNode    = hdf5->groups;
4455c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
4465c6c1daeSBarry Smith   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
4475c6c1daeSBarry Smith   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
4485c6c1daeSBarry Smith   PetscFunctionReturn(0);
4495c6c1daeSBarry Smith }
4505c6c1daeSBarry Smith 
4515c6c1daeSBarry Smith #undef __FUNCT__
4525c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5GetGroup"
4535c6c1daeSBarry Smith /*@C
4540298fd71SBarry Smith   PetscViewerHDF5GetGroup - Get the current HDF5 group for output. If none has been assigned, returns NULL.
4555c6c1daeSBarry Smith 
4565c6c1daeSBarry Smith   Not collective
4575c6c1daeSBarry Smith 
4585c6c1daeSBarry Smith   Input Parameter:
4595c6c1daeSBarry Smith . viewer - the PetscViewer
4605c6c1daeSBarry Smith 
4615c6c1daeSBarry Smith   Output Parameter:
4625c6c1daeSBarry Smith . name - The group name
4635c6c1daeSBarry Smith 
4645c6c1daeSBarry Smith   Level: intermediate
4655c6c1daeSBarry Smith 
4665c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup()
4675c6c1daeSBarry Smith @*/
4685c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
4695c6c1daeSBarry Smith {
4705c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
4715c6c1daeSBarry Smith 
4725c6c1daeSBarry Smith   PetscFunctionBegin;
4735c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
474c959eef4SJed Brown   PetscValidPointer(name,2);
475a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
4760298fd71SBarry Smith   else *name = NULL;
4775c6c1daeSBarry Smith   PetscFunctionReturn(0);
4785c6c1daeSBarry Smith }
4795c6c1daeSBarry Smith 
4805c6c1daeSBarry Smith #undef __FUNCT__
4815c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5IncrementTimestep"
4823ef9c667SSatish Balay /*@
4835c6c1daeSBarry Smith   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
4845c6c1daeSBarry Smith 
4855c6c1daeSBarry Smith   Not collective
4865c6c1daeSBarry Smith 
4875c6c1daeSBarry Smith   Input Parameter:
4885c6c1daeSBarry Smith . viewer - the PetscViewer
4895c6c1daeSBarry Smith 
4905c6c1daeSBarry Smith   Level: intermediate
4915c6c1daeSBarry Smith 
4925c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
4935c6c1daeSBarry Smith @*/
4945c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
4955c6c1daeSBarry Smith {
4965c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4975c6c1daeSBarry Smith 
4985c6c1daeSBarry Smith   PetscFunctionBegin;
4995c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5005c6c1daeSBarry Smith   ++hdf5->timestep;
5015c6c1daeSBarry Smith   PetscFunctionReturn(0);
5025c6c1daeSBarry Smith }
5035c6c1daeSBarry Smith 
5045c6c1daeSBarry Smith #undef __FUNCT__
5055c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5SetTimestep"
5063ef9c667SSatish Balay /*@
5075c6c1daeSBarry Smith   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
5085c6c1daeSBarry Smith   of -1 disables blocking with timesteps.
5095c6c1daeSBarry Smith 
5105c6c1daeSBarry Smith   Not collective
5115c6c1daeSBarry Smith 
5125c6c1daeSBarry Smith   Input Parameters:
5135c6c1daeSBarry Smith + viewer - the PetscViewer
5145c6c1daeSBarry Smith - timestep - The timestep number
5155c6c1daeSBarry Smith 
5165c6c1daeSBarry Smith   Level: intermediate
5175c6c1daeSBarry Smith 
5185c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
5195c6c1daeSBarry Smith @*/
5205c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
5215c6c1daeSBarry Smith {
5225c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5235c6c1daeSBarry Smith 
5245c6c1daeSBarry Smith   PetscFunctionBegin;
5255c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5265c6c1daeSBarry Smith   hdf5->timestep = timestep;
5275c6c1daeSBarry Smith   PetscFunctionReturn(0);
5285c6c1daeSBarry Smith }
5295c6c1daeSBarry Smith 
5305c6c1daeSBarry Smith #undef __FUNCT__
5315c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5GetTimestep"
5323ef9c667SSatish Balay /*@
5335c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
5345c6c1daeSBarry Smith 
5355c6c1daeSBarry Smith   Not collective
5365c6c1daeSBarry Smith 
5375c6c1daeSBarry Smith   Input Parameter:
5385c6c1daeSBarry Smith . viewer - the PetscViewer
5395c6c1daeSBarry Smith 
5405c6c1daeSBarry Smith   Output Parameter:
5415c6c1daeSBarry Smith . timestep - The timestep number
5425c6c1daeSBarry Smith 
5435c6c1daeSBarry Smith   Level: intermediate
5445c6c1daeSBarry Smith 
5455c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
5465c6c1daeSBarry Smith @*/
5475c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
5485c6c1daeSBarry Smith {
5495c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5505c6c1daeSBarry Smith 
5515c6c1daeSBarry Smith   PetscFunctionBegin;
5525c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5535c6c1daeSBarry Smith   PetscValidPointer(timestep,2);
5545c6c1daeSBarry Smith   *timestep = hdf5->timestep;
5555c6c1daeSBarry Smith   PetscFunctionReturn(0);
5565c6c1daeSBarry Smith }
5575c6c1daeSBarry Smith 
55836bce6e8SMatthew G. Knepley #undef __FUNCT__
55936bce6e8SMatthew G. Knepley #define __FUNCT__ "PetscDataTypeToHDF5DataType"
56036bce6e8SMatthew G. Knepley /*@C
56136bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
56236bce6e8SMatthew G. Knepley 
56336bce6e8SMatthew G. Knepley   Not collective
56436bce6e8SMatthew G. Knepley 
56536bce6e8SMatthew G. Knepley   Input Parameter:
56636bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
56736bce6e8SMatthew G. Knepley 
56836bce6e8SMatthew G. Knepley   Output Parameter:
56936bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
57036bce6e8SMatthew G. Knepley 
57136bce6e8SMatthew G. Knepley   Level: advanced
57236bce6e8SMatthew G. Knepley 
57336bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
57436bce6e8SMatthew G. Knepley @*/
57536bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
57636bce6e8SMatthew G. Knepley {
57736bce6e8SMatthew G. Knepley   PetscFunctionBegin;
57836bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
57936bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
58036bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_LLONG;
58136bce6e8SMatthew G. Knepley #else
58236bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_INT;
58336bce6e8SMatthew G. Knepley #endif
58436bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
58536bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
58636bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
58736bce6e8SMatthew G. Knepley   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_DOUBLE;
58836bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_DOUBLE;
58936bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
59036bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
59136bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
5927e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
59336bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
59436bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
59536bce6e8SMatthew G. Knepley }
59636bce6e8SMatthew G. Knepley 
59736bce6e8SMatthew G. Knepley #undef __FUNCT__
59836bce6e8SMatthew G. Knepley #define __FUNCT__ "PetscHDF5DataTypeToPetscDataType"
59936bce6e8SMatthew G. Knepley /*@C
60036bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
60136bce6e8SMatthew G. Knepley 
60236bce6e8SMatthew G. Knepley   Not collective
60336bce6e8SMatthew G. Knepley 
60436bce6e8SMatthew G. Knepley   Input Parameter:
60536bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
60636bce6e8SMatthew G. Knepley 
60736bce6e8SMatthew G. Knepley   Output Parameter:
60836bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
60936bce6e8SMatthew G. Knepley 
61036bce6e8SMatthew G. Knepley   Level: advanced
61136bce6e8SMatthew G. Knepley 
61236bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
61336bce6e8SMatthew G. Knepley @*/
61436bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
61536bce6e8SMatthew G. Knepley {
61636bce6e8SMatthew G. Knepley   PetscFunctionBegin;
61736bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
61836bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
61936bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
62036bce6e8SMatthew G. Knepley #else
62136bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
62236bce6e8SMatthew G. Knepley #endif
62336bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
62436bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
62536bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
62636bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
62736bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
62836bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
6297e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
63036bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
63136bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
63236bce6e8SMatthew G. Knepley }
63336bce6e8SMatthew G. Knepley 
63436bce6e8SMatthew G. Knepley #undef __FUNCT__
63536bce6e8SMatthew G. Knepley #define __FUNCT__ "PetscViewerHDF5WriteAttribute"
636df863907SAlex Fikl /*@C
63736bce6e8SMatthew G. Knepley  PetscViewerHDF5WriteAttribute - Write a scalar attribute
63836bce6e8SMatthew G. Knepley 
63936bce6e8SMatthew G. Knepley   Input Parameters:
64036bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer
64136bce6e8SMatthew G. Knepley . parent - The parent name
64236bce6e8SMatthew G. Knepley . name   - The attribute name
64336bce6e8SMatthew G. Knepley . datatype - The attribute type
64436bce6e8SMatthew G. Knepley - value    - The attribute value
64536bce6e8SMatthew G. Knepley 
64636bce6e8SMatthew G. Knepley   Level: advanced
64736bce6e8SMatthew G. Knepley 
648e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute()
64936bce6e8SMatthew G. Knepley @*/
65036bce6e8SMatthew G. Knepley PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
65136bce6e8SMatthew G. Knepley {
652729ab6d9SBarry Smith   hid_t          h5, dataspace, dataset, attribute, dtype;
65336bce6e8SMatthew G. Knepley   PetscErrorCode ierr;
65436bce6e8SMatthew G. Knepley 
65536bce6e8SMatthew G. Knepley   PetscFunctionBegin;
6565cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
65736bce6e8SMatthew G. Knepley   PetscValidPointer(parent, 2);
65836bce6e8SMatthew G. Knepley   PetscValidPointer(name, 3);
65936bce6e8SMatthew G. Knepley   PetscValidPointer(value, 4);
66036bce6e8SMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
6617e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
6627e97c476SMatthew G. Knepley     size_t len;
6637e97c476SMatthew G. Knepley     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
664729ab6d9SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
6657e97c476SMatthew G. Knepley   }
66636bce6e8SMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
667729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
66836bce6e8SMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
669729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataset,H5Dopen2,(h5, parent, H5P_DEFAULT));
670729ab6d9SBarry Smith   PetscStackCallHDF5Return(attribute,H5Acreate2,(dataset, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
67136bce6e8SMatthew G. Knepley #else
672729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataset,H5Dopen,(h5, parent));
673729ab6d9SBarry Smith   PetscStackCallHDF5Return(attribute,H5Acreate,(dataset, name, dtype, dataspace, H5P_DEFAULT));
67436bce6e8SMatthew G. Knepley #endif
675729ab6d9SBarry Smith   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
676729ab6d9SBarry Smith   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
677729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
678729ab6d9SBarry Smith   PetscStackCallHDF5(H5Dclose,(dataset));
679729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(dataspace));
68036bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
68136bce6e8SMatthew G. Knepley }
68236bce6e8SMatthew G. Knepley 
683ad0c61feSMatthew G. Knepley #undef __FUNCT__
684ad0c61feSMatthew G. Knepley #define __FUNCT__ "PetscViewerHDF5ReadAttribute"
685df863907SAlex Fikl /*@C
686ad0c61feSMatthew G. Knepley  PetscViewerHDF5ReadAttribute - Read a scalar attribute
687ad0c61feSMatthew G. Knepley 
688ad0c61feSMatthew G. Knepley   Input Parameters:
689ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer
690ad0c61feSMatthew G. Knepley . parent - The parent name
691ad0c61feSMatthew G. Knepley . name   - The attribute name
692ad0c61feSMatthew G. Knepley - datatype - The attribute type
693ad0c61feSMatthew G. Knepley 
694ad0c61feSMatthew G. Knepley   Output Parameter:
695ad0c61feSMatthew G. Knepley . value    - The attribute value
696ad0c61feSMatthew G. Knepley 
697ad0c61feSMatthew G. Knepley   Level: advanced
698ad0c61feSMatthew G. Knepley 
699e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute()
700ad0c61feSMatthew G. Knepley @*/
701ad0c61feSMatthew G. Knepley PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, void *value)
702ad0c61feSMatthew G. Knepley {
703729ab6d9SBarry Smith   hid_t          h5, dataspace, dataset, attribute, dtype;
704ad0c61feSMatthew G. Knepley   PetscErrorCode ierr;
705ad0c61feSMatthew G. Knepley 
706ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
7075cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
708ad0c61feSMatthew G. Knepley   PetscValidPointer(parent, 2);
709ad0c61feSMatthew G. Knepley   PetscValidPointer(name, 3);
710ad0c61feSMatthew G. Knepley   PetscValidPointer(value, 4);
711ad0c61feSMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
712ad0c61feSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
713729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
714ad0c61feSMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
715729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataset,H5Dopen2,(h5, parent, H5P_DEFAULT));
716ad0c61feSMatthew G. Knepley #else
717729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataset,H5Dopen,(h5, parent));
718ad0c61feSMatthew G. Knepley #endif
719729ab6d9SBarry Smith   PetscStackCallHDF5Return(attribute,H5Aopen_name,(dataset, name));
72070efba33SBarry Smith   PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
721729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
722729ab6d9SBarry Smith   PetscStackCallHDF5(H5Dclose,(dataset));
723729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(dataspace));
724ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
725ad0c61feSMatthew G. Knepley }
726ad0c61feSMatthew G. Knepley 
727e7bdbf8eSMatthew G. Knepley #undef __FUNCT__
7285cdeb108SMatthew G. Knepley #define __FUNCT__ "PetscViewerHDF5HasObject"
7295cdeb108SMatthew G. Knepley static PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, const char name[], H5O_type_t otype, PetscBool *has)
7305cdeb108SMatthew G. Knepley {
7315cdeb108SMatthew G. Knepley   hid_t          h5;
7325cdeb108SMatthew G. Knepley   PetscErrorCode ierr;
7335cdeb108SMatthew G. Knepley 
7345cdeb108SMatthew G. Knepley   PetscFunctionBegin;
7355cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7365cdeb108SMatthew G. Knepley   PetscValidPointer(name, 2);
7375cdeb108SMatthew G. Knepley   PetscValidPointer(has, 3);
7385cdeb108SMatthew G. Knepley   *has = PETSC_FALSE;
739c0cd0301SJed Brown   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
7405cdeb108SMatthew G. Knepley   if (H5Lexists(h5, name, H5P_DEFAULT)) {
7415cdeb108SMatthew G. Knepley     H5O_info_t info;
7425cdeb108SMatthew G. Knepley     hid_t      obj;
7435cdeb108SMatthew G. Knepley 
744729ab6d9SBarry Smith     PetscStackCallHDF5Return(obj,H5Oopen,(h5, name, H5P_DEFAULT));
745729ab6d9SBarry Smith     PetscStackCallHDF5(H5Oget_info,(obj, &info));
7465cdeb108SMatthew G. Knepley     if (otype == info.type) *has = PETSC_TRUE;
747729ab6d9SBarry Smith     PetscStackCallHDF5(H5Oclose,(obj));
7485cdeb108SMatthew G. Knepley   }
7495cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
7505cdeb108SMatthew G. Knepley }
7515cdeb108SMatthew G. Knepley 
7525cdeb108SMatthew G. Knepley #undef __FUNCT__
753e7bdbf8eSMatthew G. Knepley #define __FUNCT__ "PetscViewerHDF5HasAttribute"
754df863907SAlex Fikl /*@C
755e7bdbf8eSMatthew G. Knepley  PetscViewerHDF5HasAttribute - Check whether a scalar attribute exists
756e7bdbf8eSMatthew G. Knepley 
757e7bdbf8eSMatthew G. Knepley   Input Parameters:
758e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer
759e7bdbf8eSMatthew G. Knepley . parent - The parent name
760e7bdbf8eSMatthew G. Knepley - name   - The attribute name
761e7bdbf8eSMatthew G. Knepley 
762e7bdbf8eSMatthew G. Knepley   Output Parameter:
763e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
764e7bdbf8eSMatthew G. Knepley 
765e7bdbf8eSMatthew G. Knepley   Level: advanced
766e7bdbf8eSMatthew G. Knepley 
767e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute()
768e7bdbf8eSMatthew G. Knepley @*/
769e7bdbf8eSMatthew G. Knepley PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
770e7bdbf8eSMatthew G. Knepley {
771729ab6d9SBarry Smith   hid_t          h5, dataset;
7725cdeb108SMatthew G. Knepley   htri_t         hhas;
7735cdeb108SMatthew G. Knepley   PetscBool      exists;
774e7bdbf8eSMatthew G. Knepley   PetscErrorCode ierr;
775e7bdbf8eSMatthew G. Knepley 
776e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
7775cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
778e7bdbf8eSMatthew G. Knepley   PetscValidPointer(parent, 2);
779e7bdbf8eSMatthew G. Knepley   PetscValidPointer(name, 3);
780e7bdbf8eSMatthew G. Knepley   PetscValidPointer(has, 4);
781e7bdbf8eSMatthew G. Knepley   *has = PETSC_FALSE;
782e7bdbf8eSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
7835cdeb108SMatthew G. Knepley   ierr = PetscViewerHDF5HasObject(viewer, parent, H5O_TYPE_DATASET, &exists);CHKERRQ(ierr);
7845cdeb108SMatthew G. Knepley   if (exists) {
785e7bdbf8eSMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
786729ab6d9SBarry Smith     PetscStackCall("H5Dopen2",dataset = H5Dopen2(h5, parent, H5P_DEFAULT));
787e7bdbf8eSMatthew G. Knepley #else
788729ab6d9SBarry Smith     PetscStackCall("H5Dopen",dataset = H5Dopen(h5, parent));
789e7bdbf8eSMatthew G. Knepley #endif
790729ab6d9SBarry Smith     if (dataset < 0) PetscFunctionReturn(0);
791729ab6d9SBarry Smith     PetscStackCall("H5Aexists",hhas = H5Aexists(dataset, name));
792729ab6d9SBarry Smith     if (hhas < 0) {
793729ab6d9SBarry Smith       PetscStackCallHDF5(H5Dclose,(dataset));
794729ab6d9SBarry Smith       PetscFunctionReturn(0);
795729ab6d9SBarry Smith     }
796729ab6d9SBarry Smith     PetscStackCallHDF5(H5Dclose,(dataset));
7975cdeb108SMatthew G. Knepley     *has = hhas ? PETSC_TRUE : PETSC_FALSE;
7985cdeb108SMatthew G. Knepley   }
799e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
800e7bdbf8eSMatthew G. Knepley }
801e7bdbf8eSMatthew G. Knepley 
802a75e6a4aSMatthew G. Knepley /*
803a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
804a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
805a75e6a4aSMatthew G. Knepley */
806a75e6a4aSMatthew G. Knepley static int Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
807a75e6a4aSMatthew G. Knepley 
808a75e6a4aSMatthew G. Knepley #undef __FUNCT__
809a75e6a4aSMatthew G. Knepley #define __FUNCT__ "PETSC_VIEWER_HDF5_"
810a75e6a4aSMatthew G. Knepley /*@C
811a75e6a4aSMatthew G. Knepley   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
812a75e6a4aSMatthew G. Knepley 
813a75e6a4aSMatthew G. Knepley   Collective on MPI_Comm
814a75e6a4aSMatthew G. Knepley 
815a75e6a4aSMatthew G. Knepley   Input Parameter:
816a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer
817a75e6a4aSMatthew G. Knepley 
818a75e6a4aSMatthew G. Knepley   Level: intermediate
819a75e6a4aSMatthew G. Knepley 
820a75e6a4aSMatthew G. Knepley   Options Database Keys:
821a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name>
822a75e6a4aSMatthew G. Knepley 
823a75e6a4aSMatthew G. Knepley   Environmental variables:
824a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME
825a75e6a4aSMatthew G. Knepley 
826a75e6a4aSMatthew G. Knepley   Notes:
827a75e6a4aSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
828a75e6a4aSMatthew G. Knepley   an error code.  The HDF5 PetscViewer is usually used in the form
829a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
830a75e6a4aSMatthew G. Knepley 
831a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
832a75e6a4aSMatthew G. Knepley @*/
833a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
834a75e6a4aSMatthew G. Knepley {
835a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
836a75e6a4aSMatthew G. Knepley   PetscBool      flg;
837a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
838a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
839a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
840a75e6a4aSMatthew G. Knepley 
841a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
842a75e6a4aSMatthew 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);}
843a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
844a75e6a4aSMatthew G. Knepley     ierr = MPI_Keyval_create(MPI_NULL_COPY_FN,MPI_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
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 = MPI_Attr_get(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
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   if (!flg) { /* PetscViewer not yet created */
850a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
851a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
852a75e6a4aSMatthew G. Knepley     if (!flg) {
853a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname,"output.h5");
854a75e6a4aSMatthew G. Knepley       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
855a75e6a4aSMatthew G. Knepley     }
856a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
857a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
858a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
859a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
860a75e6a4aSMatthew G. Knepley     ierr = MPI_Attr_put(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
861a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
862a75e6a4aSMatthew G. Knepley   }
863a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
864a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
865a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
866a75e6a4aSMatthew G. Knepley }
867