xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 54dbf7064d48b386dc5817276e211a0dc72ba371)
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 
194416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v)
2082ea9b62SBarry Smith {
2182ea9b62SBarry Smith   PetscErrorCode   ierr;
2282ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;
2382ea9b62SBarry Smith 
2482ea9b62SBarry Smith   PetscFunctionBegin;
2582ea9b62SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr);
2682ea9b62SBarry Smith   ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr);
279a0502c6SHåkon Strandenes   ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr);
2882ea9b62SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
2982ea9b62SBarry Smith   PetscFunctionReturn(0);
3082ea9b62SBarry Smith }
3182ea9b62SBarry Smith 
325c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
335c6c1daeSBarry Smith {
345c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
355c6c1daeSBarry Smith   PetscErrorCode   ierr;
365c6c1daeSBarry Smith 
375c6c1daeSBarry Smith   PetscFunctionBegin;
385c6c1daeSBarry Smith   ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);
39729ab6d9SBarry Smith   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
405c6c1daeSBarry Smith   PetscFunctionReturn(0);
415c6c1daeSBarry Smith }
425c6c1daeSBarry Smith 
435c6c1daeSBarry Smith PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
445c6c1daeSBarry Smith {
455c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
465c6c1daeSBarry Smith   PetscErrorCode   ierr;
475c6c1daeSBarry Smith 
485c6c1daeSBarry Smith   PetscFunctionBegin;
495c6c1daeSBarry Smith   ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr);
505c6c1daeSBarry Smith   while (hdf5->groups) {
515c6c1daeSBarry Smith     GroupList *tmp = hdf5->groups->next;
525c6c1daeSBarry Smith 
535c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups->name);CHKERRQ(ierr);
545c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups);CHKERRQ(ierr);
555c6c1daeSBarry Smith     hdf5->groups = tmp;
565c6c1daeSBarry Smith   }
575c6c1daeSBarry Smith   ierr = PetscFree(hdf5);CHKERRQ(ierr);
580b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
59d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr);
600b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
615c6c1daeSBarry Smith   PetscFunctionReturn(0);
625c6c1daeSBarry Smith }
635c6c1daeSBarry Smith 
645c6c1daeSBarry Smith PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
655c6c1daeSBarry Smith {
665c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
675c6c1daeSBarry Smith 
685c6c1daeSBarry Smith   PetscFunctionBegin;
695c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
705c6c1daeSBarry Smith   hdf5->btype = type;
715c6c1daeSBarry Smith   PetscFunctionReturn(0);
725c6c1daeSBarry Smith }
735c6c1daeSBarry Smith 
7482ea9b62SBarry Smith PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
7582ea9b62SBarry Smith {
7682ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7782ea9b62SBarry Smith 
7882ea9b62SBarry Smith   PetscFunctionBegin;
7982ea9b62SBarry Smith   hdf5->basedimension2 = flg;
8082ea9b62SBarry Smith   PetscFunctionReturn(0);
8182ea9b62SBarry Smith }
8282ea9b62SBarry Smith 
83df863907SAlex Fikl /*@
8482ea9b62SBarry Smith      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
8582ea9b62SBarry Smith        dimension of 2.
8682ea9b62SBarry Smith 
8782ea9b62SBarry Smith     Logically Collective on PetscViewer
8882ea9b62SBarry Smith 
8982ea9b62SBarry Smith   Input Parameters:
9082ea9b62SBarry Smith +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
9182ea9b62SBarry 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
9282ea9b62SBarry Smith 
9382ea9b62SBarry Smith   Options Database:
9482ea9b62SBarry 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
9582ea9b62SBarry Smith 
9682ea9b62SBarry Smith 
9795452b02SPatrick Sanan   Notes:
9895452b02SPatrick Sanan     Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
9982ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
10082ea9b62SBarry Smith 
10182ea9b62SBarry Smith   Level: intermediate
10282ea9b62SBarry Smith 
10382ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
10482ea9b62SBarry Smith 
10582ea9b62SBarry Smith @*/
10682ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
10782ea9b62SBarry Smith {
10882ea9b62SBarry Smith   PetscErrorCode ierr;
10982ea9b62SBarry Smith 
11082ea9b62SBarry Smith   PetscFunctionBegin;
11182ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
11282ea9b62SBarry Smith   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
11382ea9b62SBarry Smith   PetscFunctionReturn(0);
11482ea9b62SBarry Smith }
11582ea9b62SBarry Smith 
116df863907SAlex Fikl /*@
11782ea9b62SBarry Smith      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
11882ea9b62SBarry Smith        dimension of 2.
11982ea9b62SBarry Smith 
12082ea9b62SBarry Smith     Logically Collective on PetscViewer
12182ea9b62SBarry Smith 
12282ea9b62SBarry Smith   Input Parameter:
12382ea9b62SBarry Smith .  viewer - the PetscViewer, must be of type HDF5
12482ea9b62SBarry Smith 
12582ea9b62SBarry Smith   Output Parameter:
12682ea9b62SBarry 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
12782ea9b62SBarry Smith 
12895452b02SPatrick Sanan   Notes:
12995452b02SPatrick Sanan     Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
13082ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
13182ea9b62SBarry Smith 
13282ea9b62SBarry Smith   Level: intermediate
13382ea9b62SBarry Smith 
13482ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
13582ea9b62SBarry Smith 
13682ea9b62SBarry Smith @*/
13782ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
13882ea9b62SBarry Smith {
13982ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
14082ea9b62SBarry Smith 
14182ea9b62SBarry Smith   PetscFunctionBegin;
14282ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
14382ea9b62SBarry Smith   *flg = hdf5->basedimension2;
14482ea9b62SBarry Smith   PetscFunctionReturn(0);
14582ea9b62SBarry Smith }
14682ea9b62SBarry Smith 
1479a0502c6SHåkon Strandenes PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
1489a0502c6SHåkon Strandenes {
1499a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1509a0502c6SHåkon Strandenes 
1519a0502c6SHåkon Strandenes   PetscFunctionBegin;
1529a0502c6SHåkon Strandenes   hdf5->spoutput = flg;
1539a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
1549a0502c6SHåkon Strandenes }
1559a0502c6SHåkon Strandenes 
156df863907SAlex Fikl /*@
1579a0502c6SHåkon Strandenes      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
1589a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
1599a0502c6SHåkon Strandenes 
1609a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
1619a0502c6SHåkon Strandenes 
1629a0502c6SHåkon Strandenes   Input Parameters:
1639a0502c6SHåkon Strandenes +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
1649a0502c6SHåkon Strandenes -  flg - if PETSC_TRUE the data will be written to disk with single precision
1659a0502c6SHåkon Strandenes 
1669a0502c6SHåkon Strandenes   Options Database:
1679a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
1689a0502c6SHåkon Strandenes 
1699a0502c6SHåkon Strandenes 
17095452b02SPatrick Sanan   Notes:
17195452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
1729a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
1739a0502c6SHåkon Strandenes 
1749a0502c6SHåkon Strandenes   Level: intermediate
1759a0502c6SHåkon Strandenes 
1769a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
1779a0502c6SHåkon Strandenes           PetscReal
1789a0502c6SHåkon Strandenes 
1799a0502c6SHåkon Strandenes @*/
1809a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
1819a0502c6SHåkon Strandenes {
1829a0502c6SHåkon Strandenes   PetscErrorCode ierr;
1839a0502c6SHåkon Strandenes 
1849a0502c6SHåkon Strandenes   PetscFunctionBegin;
1859a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1869a0502c6SHåkon Strandenes   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
1879a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
1889a0502c6SHåkon Strandenes }
1899a0502c6SHåkon Strandenes 
190df863907SAlex Fikl /*@
1919a0502c6SHåkon Strandenes      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
1929a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
1939a0502c6SHåkon Strandenes 
1949a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
1959a0502c6SHåkon Strandenes 
1969a0502c6SHåkon Strandenes   Input Parameter:
1979a0502c6SHåkon Strandenes .  viewer - the PetscViewer, must be of type HDF5
1989a0502c6SHåkon Strandenes 
1999a0502c6SHåkon Strandenes   Output Parameter:
2009a0502c6SHåkon Strandenes .  flg - if PETSC_TRUE the data will be written to disk with single precision
2019a0502c6SHåkon Strandenes 
20295452b02SPatrick Sanan   Notes:
20395452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2049a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2059a0502c6SHåkon Strandenes 
2069a0502c6SHåkon Strandenes   Level: intermediate
2079a0502c6SHåkon Strandenes 
2089a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2099a0502c6SHåkon Strandenes           PetscReal
2109a0502c6SHåkon Strandenes 
2119a0502c6SHåkon Strandenes @*/
2129a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
2139a0502c6SHåkon Strandenes {
2149a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2159a0502c6SHåkon Strandenes 
2169a0502c6SHåkon Strandenes   PetscFunctionBegin;
2179a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2189a0502c6SHåkon Strandenes   *flg = hdf5->spoutput;
2199a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2209a0502c6SHåkon Strandenes }
2219a0502c6SHåkon Strandenes 
2225c6c1daeSBarry Smith PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
2235c6c1daeSBarry Smith {
2245c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2255c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
2265c6c1daeSBarry Smith   MPI_Info          info = MPI_INFO_NULL;
2275c6c1daeSBarry Smith #endif
2285c6c1daeSBarry Smith   hid_t             plist_id;
2295c6c1daeSBarry Smith   PetscErrorCode    ierr;
2305c6c1daeSBarry Smith 
2315c6c1daeSBarry Smith   PetscFunctionBegin;
232f683cc58SToby Isaac   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
233f683cc58SToby Isaac   if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);}
2345c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
2355c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
236729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
2375c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
238729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info));
2395c6c1daeSBarry Smith #endif
2405c6c1daeSBarry Smith   /* Create or open the file collectively */
2415c6c1daeSBarry Smith   switch (hdf5->btype) {
2425c6c1daeSBarry Smith   case FILE_MODE_READ:
243729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
2445c6c1daeSBarry Smith     break;
2455c6c1daeSBarry Smith   case FILE_MODE_APPEND:
246729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
2475c6c1daeSBarry Smith     break;
2485c6c1daeSBarry Smith   case FILE_MODE_WRITE:
249729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
2505c6c1daeSBarry Smith     break;
2515c6c1daeSBarry Smith   default:
2525c6c1daeSBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
2535c6c1daeSBarry Smith   }
2545c6c1daeSBarry Smith   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
255729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
2565c6c1daeSBarry Smith   PetscFunctionReturn(0);
2575c6c1daeSBarry Smith }
2585c6c1daeSBarry Smith 
259d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
260d1232d7fSToby Isaac {
261d1232d7fSToby Isaac   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
262d1232d7fSToby Isaac 
263d1232d7fSToby Isaac   PetscFunctionBegin;
264d1232d7fSToby Isaac   *name = vhdf5->filename;
265d1232d7fSToby Isaac   PetscFunctionReturn(0);
266d1232d7fSToby Isaac }
267d1232d7fSToby Isaac 
2688556b5ebSBarry Smith /*MC
2698556b5ebSBarry Smith    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
2708556b5ebSBarry Smith 
2718556b5ebSBarry Smith 
2728556b5ebSBarry Smith .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
2738556b5ebSBarry Smith            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
2748556b5ebSBarry Smith            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
2758556b5ebSBarry Smith            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
2768556b5ebSBarry Smith 
2771b266c99SBarry Smith   Level: beginner
2788556b5ebSBarry Smith M*/
279d1232d7fSToby Isaac 
2808cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
2815c6c1daeSBarry Smith {
2825c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
2835c6c1daeSBarry Smith   PetscErrorCode   ierr;
2845c6c1daeSBarry Smith 
2855c6c1daeSBarry Smith   PetscFunctionBegin;
286b00a9115SJed Brown   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
2875c6c1daeSBarry Smith 
2885c6c1daeSBarry Smith   v->data                = (void*) hdf5;
2895c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
29082ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
2915c6c1daeSBarry Smith   v->ops->flush          = 0;
2925c6c1daeSBarry Smith   hdf5->btype            = (PetscFileMode) -1;
2935c6c1daeSBarry Smith   hdf5->filename         = 0;
2945c6c1daeSBarry Smith   hdf5->timestep         = -1;
2950298fd71SBarry Smith   hdf5->groups           = NULL;
2965c6c1daeSBarry Smith 
2970b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
298d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr);
2990b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
30082ea9b62SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
3019a0502c6SHåkon Strandenes   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr);
3025c6c1daeSBarry Smith   PetscFunctionReturn(0);
3035c6c1daeSBarry Smith }
3045c6c1daeSBarry Smith 
3055c6c1daeSBarry Smith /*@C
3065c6c1daeSBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
3075c6c1daeSBarry Smith 
3085c6c1daeSBarry Smith    Collective on MPI_Comm
3095c6c1daeSBarry Smith 
3105c6c1daeSBarry Smith    Input Parameters:
3115c6c1daeSBarry Smith +  comm - MPI communicator
3125c6c1daeSBarry Smith .  name - name of file
3135c6c1daeSBarry Smith -  type - type of file
3145c6c1daeSBarry Smith $    FILE_MODE_WRITE - create new file for binary output
3155c6c1daeSBarry Smith $    FILE_MODE_READ - open existing file for binary input
3165c6c1daeSBarry Smith $    FILE_MODE_APPEND - open existing file for binary output
3175c6c1daeSBarry Smith 
3185c6c1daeSBarry Smith    Output Parameter:
3195c6c1daeSBarry Smith .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
3205c6c1daeSBarry Smith 
32182ea9b62SBarry Smith   Options Database:
32282ea9b62SBarry 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
3239a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
32482ea9b62SBarry Smith 
3255c6c1daeSBarry Smith    Level: beginner
3265c6c1daeSBarry Smith 
3275c6c1daeSBarry Smith    Note:
3285c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
3295c6c1daeSBarry Smith 
3305c6c1daeSBarry Smith    Concepts: HDF5 files
3315c6c1daeSBarry Smith    Concepts: PetscViewerHDF5^creating
3325c6c1daeSBarry Smith 
3336a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
3349a0502c6SHåkon Strandenes           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
335a56f64adSBarry Smith           MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
3365c6c1daeSBarry Smith @*/
3375c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
3385c6c1daeSBarry Smith {
3395c6c1daeSBarry Smith   PetscErrorCode ierr;
3405c6c1daeSBarry Smith 
3415c6c1daeSBarry Smith   PetscFunctionBegin;
3425c6c1daeSBarry Smith   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
3435c6c1daeSBarry Smith   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
3445c6c1daeSBarry Smith   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
3455c6c1daeSBarry Smith   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
3465c6c1daeSBarry Smith   PetscFunctionReturn(0);
3475c6c1daeSBarry Smith }
3485c6c1daeSBarry Smith 
3495c6c1daeSBarry Smith /*@C
3505c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
3515c6c1daeSBarry Smith 
3525c6c1daeSBarry Smith   Not collective
3535c6c1daeSBarry Smith 
3545c6c1daeSBarry Smith   Input Parameter:
3555c6c1daeSBarry Smith . viewer - the PetscViewer
3565c6c1daeSBarry Smith 
3575c6c1daeSBarry Smith   Output Parameter:
3585c6c1daeSBarry Smith . file_id - The file id
3595c6c1daeSBarry Smith 
3605c6c1daeSBarry Smith   Level: intermediate
3615c6c1daeSBarry Smith 
3625c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open()
3635c6c1daeSBarry Smith @*/
3645c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
3655c6c1daeSBarry Smith {
3665c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3675c6c1daeSBarry Smith 
3685c6c1daeSBarry Smith   PetscFunctionBegin;
3695c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
3705c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
3715c6c1daeSBarry Smith   PetscFunctionReturn(0);
3725c6c1daeSBarry Smith }
3735c6c1daeSBarry Smith 
3745c6c1daeSBarry Smith /*@C
3755c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
3765c6c1daeSBarry Smith 
3775c6c1daeSBarry Smith   Not collective
3785c6c1daeSBarry Smith 
3795c6c1daeSBarry Smith   Input Parameters:
3805c6c1daeSBarry Smith + viewer - the PetscViewer
3815c6c1daeSBarry Smith - name - The group name
3825c6c1daeSBarry Smith 
3835c6c1daeSBarry Smith   Level: intermediate
3845c6c1daeSBarry Smith 
3855c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
3865c6c1daeSBarry Smith @*/
3875c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
3885c6c1daeSBarry Smith {
3895c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3905c6c1daeSBarry Smith   GroupList        *groupNode;
3915c6c1daeSBarry Smith   PetscErrorCode   ierr;
3925c6c1daeSBarry Smith 
3935c6c1daeSBarry Smith   PetscFunctionBegin;
3945c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
3955c6c1daeSBarry Smith   PetscValidCharPointer(name,2);
39695dccacaSBarry Smith   ierr = PetscNew(&groupNode);CHKERRQ(ierr);
3975c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
398a297a907SKarl Rupp 
3995c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
4005c6c1daeSBarry Smith   hdf5->groups    = groupNode;
4015c6c1daeSBarry Smith   PetscFunctionReturn(0);
4025c6c1daeSBarry Smith }
4035c6c1daeSBarry Smith 
4043ef9c667SSatish Balay /*@
4055c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
4065c6c1daeSBarry Smith 
4075c6c1daeSBarry Smith   Not collective
4085c6c1daeSBarry Smith 
4095c6c1daeSBarry Smith   Input Parameter:
4105c6c1daeSBarry Smith . viewer - the PetscViewer
4115c6c1daeSBarry Smith 
4125c6c1daeSBarry Smith   Level: intermediate
4135c6c1daeSBarry Smith 
4145c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup()
4155c6c1daeSBarry Smith @*/
4165c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
4175c6c1daeSBarry Smith {
4185c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4195c6c1daeSBarry Smith   GroupList        *groupNode;
4205c6c1daeSBarry Smith   PetscErrorCode   ierr;
4215c6c1daeSBarry Smith 
4225c6c1daeSBarry Smith   PetscFunctionBegin;
4235c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
42482f516ccSBarry Smith   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
4255c6c1daeSBarry Smith   groupNode    = hdf5->groups;
4265c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
4275c6c1daeSBarry Smith   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
4285c6c1daeSBarry Smith   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
4295c6c1daeSBarry Smith   PetscFunctionReturn(0);
4305c6c1daeSBarry Smith }
4315c6c1daeSBarry Smith 
4325c6c1daeSBarry Smith /*@C
4330298fd71SBarry Smith   PetscViewerHDF5GetGroup - Get the current HDF5 group for output. If none has been assigned, returns NULL.
4345c6c1daeSBarry Smith 
4355c6c1daeSBarry Smith   Not collective
4365c6c1daeSBarry Smith 
4375c6c1daeSBarry Smith   Input Parameter:
4385c6c1daeSBarry Smith . viewer - the PetscViewer
4395c6c1daeSBarry Smith 
4405c6c1daeSBarry Smith   Output Parameter:
4415c6c1daeSBarry Smith . name - The group name
4425c6c1daeSBarry Smith 
4435c6c1daeSBarry Smith   Level: intermediate
4445c6c1daeSBarry Smith 
4455c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup()
4465c6c1daeSBarry Smith @*/
4475c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
4485c6c1daeSBarry Smith {
4495c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
4505c6c1daeSBarry Smith 
4515c6c1daeSBarry Smith   PetscFunctionBegin;
4525c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
453c959eef4SJed Brown   PetscValidPointer(name,2);
454a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
4550298fd71SBarry Smith   else *name = NULL;
4565c6c1daeSBarry Smith   PetscFunctionReturn(0);
4575c6c1daeSBarry Smith }
4585c6c1daeSBarry Smith 
459*54dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
460*54dbf706SVaclav Hapla {
461*54dbf706SVaclav Hapla   hid_t          file_id, group;
462*54dbf706SVaclav Hapla   htri_t         found;
463*54dbf706SVaclav Hapla   const char     *groupName = NULL;
464*54dbf706SVaclav Hapla   PetscErrorCode ierr;
465*54dbf706SVaclav Hapla 
466*54dbf706SVaclav Hapla   PetscFunctionBegin;
467*54dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr);
468*54dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr);
469*54dbf706SVaclav Hapla   /* Open group */
470*54dbf706SVaclav Hapla   if (groupName) {
471*54dbf706SVaclav Hapla     PetscBool root;
472*54dbf706SVaclav Hapla 
473*54dbf706SVaclav Hapla     ierr = PetscStrcmp(groupName, "/", &root);CHKERRQ(ierr);
474*54dbf706SVaclav Hapla     PetscStackCall("H5Lexists",found = H5Lexists(file_id, groupName, H5P_DEFAULT));
475*54dbf706SVaclav Hapla     if (!root && (found <= 0)) {
476*54dbf706SVaclav Hapla #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
477*54dbf706SVaclav Hapla       PetscStackCallHDF5Return(group,H5Gcreate2,(file_id, groupName, 0, H5P_DEFAULT, H5P_DEFAULT));
478*54dbf706SVaclav Hapla #else /* deprecated HDF5 1.6 API */
479*54dbf706SVaclav Hapla       PetscStackCallHDF5Return(group,H5Gcreate,(file_id, groupName, 0));
480*54dbf706SVaclav Hapla #endif
481*54dbf706SVaclav Hapla       PetscStackCallHDF5(H5Gclose,(group));
482*54dbf706SVaclav Hapla     }
483*54dbf706SVaclav Hapla #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
484*54dbf706SVaclav Hapla     PetscStackCallHDF5Return(group,H5Gopen2,(file_id, groupName, H5P_DEFAULT));
485*54dbf706SVaclav Hapla #else
486*54dbf706SVaclav Hapla     PetscStackCallHDF5Return(group,H5Gopen,(file_id, groupName));
487*54dbf706SVaclav Hapla #endif
488*54dbf706SVaclav Hapla   } else group = file_id;
489*54dbf706SVaclav Hapla 
490*54dbf706SVaclav Hapla   *fileId  = file_id;
491*54dbf706SVaclav Hapla   *groupId = group;
492*54dbf706SVaclav Hapla   PetscFunctionReturn(0);
493*54dbf706SVaclav Hapla }
494*54dbf706SVaclav Hapla 
4953ef9c667SSatish Balay /*@
4965c6c1daeSBarry Smith   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
4975c6c1daeSBarry Smith 
4985c6c1daeSBarry Smith   Not collective
4995c6c1daeSBarry Smith 
5005c6c1daeSBarry Smith   Input Parameter:
5015c6c1daeSBarry Smith . viewer - the PetscViewer
5025c6c1daeSBarry Smith 
5035c6c1daeSBarry Smith   Level: intermediate
5045c6c1daeSBarry Smith 
5055c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
5065c6c1daeSBarry Smith @*/
5075c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
5085c6c1daeSBarry Smith {
5095c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5105c6c1daeSBarry Smith 
5115c6c1daeSBarry Smith   PetscFunctionBegin;
5125c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5135c6c1daeSBarry Smith   ++hdf5->timestep;
5145c6c1daeSBarry Smith   PetscFunctionReturn(0);
5155c6c1daeSBarry Smith }
5165c6c1daeSBarry Smith 
5173ef9c667SSatish Balay /*@
5185c6c1daeSBarry Smith   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
5195c6c1daeSBarry Smith   of -1 disables blocking with timesteps.
5205c6c1daeSBarry Smith 
5215c6c1daeSBarry Smith   Not collective
5225c6c1daeSBarry Smith 
5235c6c1daeSBarry Smith   Input Parameters:
5245c6c1daeSBarry Smith + viewer - the PetscViewer
5255c6c1daeSBarry Smith - timestep - The timestep number
5265c6c1daeSBarry Smith 
5275c6c1daeSBarry Smith   Level: intermediate
5285c6c1daeSBarry Smith 
5295c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
5305c6c1daeSBarry Smith @*/
5315c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5SetTimestep(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   hdf5->timestep = timestep;
5385c6c1daeSBarry Smith   PetscFunctionReturn(0);
5395c6c1daeSBarry Smith }
5405c6c1daeSBarry Smith 
5413ef9c667SSatish Balay /*@
5425c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
5435c6c1daeSBarry Smith 
5445c6c1daeSBarry Smith   Not collective
5455c6c1daeSBarry Smith 
5465c6c1daeSBarry Smith   Input Parameter:
5475c6c1daeSBarry Smith . viewer - the PetscViewer
5485c6c1daeSBarry Smith 
5495c6c1daeSBarry Smith   Output Parameter:
5505c6c1daeSBarry Smith . timestep - The timestep number
5515c6c1daeSBarry Smith 
5525c6c1daeSBarry Smith   Level: intermediate
5535c6c1daeSBarry Smith 
5545c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
5555c6c1daeSBarry Smith @*/
5565c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
5575c6c1daeSBarry Smith {
5585c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5595c6c1daeSBarry Smith 
5605c6c1daeSBarry Smith   PetscFunctionBegin;
5615c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5625c6c1daeSBarry Smith   PetscValidPointer(timestep,2);
5635c6c1daeSBarry Smith   *timestep = hdf5->timestep;
5645c6c1daeSBarry Smith   PetscFunctionReturn(0);
5655c6c1daeSBarry Smith }
5665c6c1daeSBarry Smith 
56736bce6e8SMatthew G. Knepley /*@C
56836bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
56936bce6e8SMatthew G. Knepley 
57036bce6e8SMatthew G. Knepley   Not collective
57136bce6e8SMatthew G. Knepley 
57236bce6e8SMatthew G. Knepley   Input Parameter:
57336bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
57436bce6e8SMatthew G. Knepley 
57536bce6e8SMatthew G. Knepley   Output Parameter:
57636bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
57736bce6e8SMatthew G. Knepley 
57836bce6e8SMatthew G. Knepley   Level: advanced
57936bce6e8SMatthew G. Knepley 
58036bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
58136bce6e8SMatthew G. Knepley @*/
58236bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
58336bce6e8SMatthew G. Knepley {
58436bce6e8SMatthew G. Knepley   PetscFunctionBegin;
58536bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
58636bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
58736bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_LLONG;
58836bce6e8SMatthew G. Knepley #else
58936bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_INT;
59036bce6e8SMatthew G. Knepley #endif
59136bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
59236bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
59336bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
59436bce6e8SMatthew G. Knepley   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_DOUBLE;
59536bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_DOUBLE;
59636bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
59736bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
59836bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
5997e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
60036bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
60136bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
60236bce6e8SMatthew G. Knepley }
60336bce6e8SMatthew G. Knepley 
60436bce6e8SMatthew G. Knepley /*@C
60536bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
60636bce6e8SMatthew G. Knepley 
60736bce6e8SMatthew G. Knepley   Not collective
60836bce6e8SMatthew G. Knepley 
60936bce6e8SMatthew G. Knepley   Input Parameter:
61036bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
61136bce6e8SMatthew G. Knepley 
61236bce6e8SMatthew G. Knepley   Output Parameter:
61336bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
61436bce6e8SMatthew G. Knepley 
61536bce6e8SMatthew G. Knepley   Level: advanced
61636bce6e8SMatthew G. Knepley 
61736bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
61836bce6e8SMatthew G. Knepley @*/
61936bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
62036bce6e8SMatthew G. Knepley {
62136bce6e8SMatthew G. Knepley   PetscFunctionBegin;
62236bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
62336bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
62436bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
62536bce6e8SMatthew G. Knepley #else
62636bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
62736bce6e8SMatthew G. Knepley #endif
62836bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
62936bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
63036bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
63136bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
63236bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
63336bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
6347e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
63536bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
63636bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
63736bce6e8SMatthew G. Knepley }
63836bce6e8SMatthew G. Knepley 
639df863907SAlex Fikl /*@C
64036bce6e8SMatthew G. Knepley  PetscViewerHDF5WriteAttribute - Write a scalar attribute
64136bce6e8SMatthew G. Knepley 
64236bce6e8SMatthew G. Knepley   Input Parameters:
64336bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer
64436bce6e8SMatthew G. Knepley . parent - The parent name
64536bce6e8SMatthew G. Knepley . name   - The attribute name
64636bce6e8SMatthew G. Knepley . datatype - The attribute type
64736bce6e8SMatthew G. Knepley - value    - The attribute value
64836bce6e8SMatthew G. Knepley 
64936bce6e8SMatthew G. Knepley   Level: advanced
65036bce6e8SMatthew G. Knepley 
651e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute()
65236bce6e8SMatthew G. Knepley @*/
65336bce6e8SMatthew G. Knepley PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
65436bce6e8SMatthew G. Knepley {
65560568592SMatthew G. Knepley   hid_t          h5, dataspace, obj, attribute, dtype;
65636bce6e8SMatthew G. Knepley   PetscErrorCode ierr;
65736bce6e8SMatthew G. Knepley 
65836bce6e8SMatthew G. Knepley   PetscFunctionBegin;
6595cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
66036bce6e8SMatthew G. Knepley   PetscValidPointer(parent, 2);
66136bce6e8SMatthew G. Knepley   PetscValidPointer(name, 3);
66236bce6e8SMatthew G. Knepley   PetscValidPointer(value, 4);
66336bce6e8SMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
6647e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
6657e97c476SMatthew G. Knepley     size_t len;
6667e97c476SMatthew G. Knepley     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
667729ab6d9SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
6687e97c476SMatthew G. Knepley   }
66936bce6e8SMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
670729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
67160568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
67236bce6e8SMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
67360568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
67436bce6e8SMatthew G. Knepley #else
67560568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Acreate,(obj, name, dtype, dataspace, H5P_DEFAULT));
67636bce6e8SMatthew G. Knepley #endif
677729ab6d9SBarry Smith   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
678729ab6d9SBarry Smith   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
679729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
68060568592SMatthew G. Knepley   PetscStackCallHDF5(H5Oclose,(obj));
681729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(dataspace));
68236bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
68336bce6e8SMatthew G. Knepley }
68436bce6e8SMatthew G. Knepley 
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 {
703f0b58479SMatthew G. Knepley   hid_t          h5, obj, attribute, atype, 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);
71360568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
71460568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
715f0b58479SMatthew G. Knepley   PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
716f0b58479SMatthew G. Knepley   if (datatype == PETSC_STRING) {
717f0b58479SMatthew G. Knepley     size_t len;
718f0b58479SMatthew G. Knepley 
719f0b58479SMatthew G. Knepley     PetscStackCallHDF5Return(len,H5Tget_size,(atype));
720f0b58479SMatthew G. Knepley     PetscStackCallHDF5(H5Tclose,(atype));
721f0b58479SMatthew G. Knepley     ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr);
722f0b58479SMatthew G. Knepley   }
72370efba33SBarry Smith   PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
724729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
72560568592SMatthew G. Knepley   PetscStackCallHDF5(H5Dclose,(obj));
726ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
727ad0c61feSMatthew G. Knepley }
728ad0c61feSMatthew G. Knepley 
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 
752df863907SAlex Fikl /*@C
753e7bdbf8eSMatthew G. Knepley  PetscViewerHDF5HasAttribute - Check whether a scalar attribute exists
754e7bdbf8eSMatthew G. Knepley 
755e7bdbf8eSMatthew G. Knepley   Input Parameters:
756e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer
757e7bdbf8eSMatthew G. Knepley . parent - The parent name
758e7bdbf8eSMatthew G. Knepley - name   - The attribute name
759e7bdbf8eSMatthew G. Knepley 
760e7bdbf8eSMatthew G. Knepley   Output Parameter:
761e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
762e7bdbf8eSMatthew G. Knepley 
763e7bdbf8eSMatthew G. Knepley   Level: advanced
764e7bdbf8eSMatthew G. Knepley 
765e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute()
766e7bdbf8eSMatthew G. Knepley @*/
767e7bdbf8eSMatthew G. Knepley PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
768e7bdbf8eSMatthew G. Knepley {
769729ab6d9SBarry Smith   hid_t          h5, dataset;
7705cdeb108SMatthew G. Knepley   htri_t         hhas;
7715cdeb108SMatthew G. Knepley   PetscBool      exists;
772e7bdbf8eSMatthew G. Knepley   PetscErrorCode ierr;
773e7bdbf8eSMatthew G. Knepley 
774e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
7755cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
776e7bdbf8eSMatthew G. Knepley   PetscValidPointer(parent, 2);
777e7bdbf8eSMatthew G. Knepley   PetscValidPointer(name, 3);
778e7bdbf8eSMatthew G. Knepley   PetscValidPointer(has, 4);
779e7bdbf8eSMatthew G. Knepley   *has = PETSC_FALSE;
780e7bdbf8eSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
7815cdeb108SMatthew G. Knepley   ierr = PetscViewerHDF5HasObject(viewer, parent, H5O_TYPE_DATASET, &exists);CHKERRQ(ierr);
7825cdeb108SMatthew G. Knepley   if (exists) {
783e7bdbf8eSMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
784729ab6d9SBarry Smith     PetscStackCall("H5Dopen2",dataset = H5Dopen2(h5, parent, H5P_DEFAULT));
785e7bdbf8eSMatthew G. Knepley #else
786729ab6d9SBarry Smith     PetscStackCall("H5Dopen",dataset = H5Dopen(h5, parent));
787e7bdbf8eSMatthew G. Knepley #endif
788729ab6d9SBarry Smith     if (dataset < 0) PetscFunctionReturn(0);
789729ab6d9SBarry Smith     PetscStackCall("H5Aexists",hhas = H5Aexists(dataset, name));
790729ab6d9SBarry Smith     if (hhas < 0) {
791729ab6d9SBarry Smith       PetscStackCallHDF5(H5Dclose,(dataset));
792729ab6d9SBarry Smith       PetscFunctionReturn(0);
793729ab6d9SBarry Smith     }
794729ab6d9SBarry Smith     PetscStackCallHDF5(H5Dclose,(dataset));
7955cdeb108SMatthew G. Knepley     *has = hhas ? PETSC_TRUE : PETSC_FALSE;
7965cdeb108SMatthew G. Knepley   }
797e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
798e7bdbf8eSMatthew G. Knepley }
799e7bdbf8eSMatthew G. Knepley 
800a75e6a4aSMatthew G. Knepley /*
801a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
802a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
803a75e6a4aSMatthew G. Knepley */
804d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
805a75e6a4aSMatthew G. Knepley 
806a75e6a4aSMatthew G. Knepley /*@C
807a75e6a4aSMatthew G. Knepley   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
808a75e6a4aSMatthew G. Knepley 
809a75e6a4aSMatthew G. Knepley   Collective on MPI_Comm
810a75e6a4aSMatthew G. Knepley 
811a75e6a4aSMatthew G. Knepley   Input Parameter:
812a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer
813a75e6a4aSMatthew G. Knepley 
814a75e6a4aSMatthew G. Knepley   Level: intermediate
815a75e6a4aSMatthew G. Knepley 
816a75e6a4aSMatthew G. Knepley   Options Database Keys:
817a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name>
818a75e6a4aSMatthew G. Knepley 
819a75e6a4aSMatthew G. Knepley   Environmental variables:
820a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME
821a75e6a4aSMatthew G. Knepley 
822a75e6a4aSMatthew G. Knepley   Notes:
823a75e6a4aSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
824a75e6a4aSMatthew G. Knepley   an error code.  The HDF5 PetscViewer is usually used in the form
825a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
826a75e6a4aSMatthew G. Knepley 
827a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
828a75e6a4aSMatthew G. Knepley @*/
829a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
830a75e6a4aSMatthew G. Knepley {
831a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
832a75e6a4aSMatthew G. Knepley   PetscBool      flg;
833a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
834a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
835a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
836a75e6a4aSMatthew G. Knepley 
837a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
838a75e6a4aSMatthew 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);}
839a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
84012801b39SBarry Smith     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
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   }
84347435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
844a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
845a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
846a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
847a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
848a75e6a4aSMatthew G. Knepley     if (!flg) {
849a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname,"output.h5");
850a75e6a4aSMatthew G. Knepley       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
851a75e6a4aSMatthew G. Knepley     }
852a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
853a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
854a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
855a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
85647435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)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   }
859a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
860a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
861a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
862a75e6a4aSMatthew G. Knepley }
863