xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision cd2571c5d4122701d9c8a8814fb939e6fb00b2ad)
1 #include <petsc/private/viewerimpl.h>
2 #include <petsc/private/viewerhdf5impl.h>
3 #include <petscviewerhdf5.h>    /*I   "petscviewerhdf5.h"   I*/
4 
5 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool*, H5O_type_t*);
6 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool*);
7 
8 static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char path[], char *abspath[])
9 {
10   PetscBool      relative = PETSC_FALSE;
11   const char     *group;
12   char           buf[PETSC_MAX_PATH_LEN] = "";
13   PetscErrorCode ierr;
14 
15   PetscFunctionBegin;
16   ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr);
17   ierr = PetscViewerHDF5PathIsRelative(path, PETSC_TRUE, &relative);CHKERRQ(ierr);
18   if (relative) {
19     ierr = PetscStrcpy(buf, group);CHKERRQ(ierr);
20     ierr = PetscStrcat(buf, "/");CHKERRQ(ierr);
21     ierr = PetscStrcat(buf, path);CHKERRQ(ierr);
22     ierr = PetscStrallocpy(buf, abspath);CHKERRQ(ierr);
23   } else {
24     ierr = PetscStrallocpy(path, abspath);CHKERRQ(ierr);
25   }
26   PetscFunctionReturn(0);
27 }
28 
29 static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj)
30 {
31   PetscBool      has;
32   PetscErrorCode ierr;
33 
34   PetscFunctionBegin;
35   ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr);
36   if (!has) {
37     const char *group;
38     ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr);
39     SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", obj->name, group ? group : "/");
40   }
41   PetscFunctionReturn(0);
42 }
43 
44 static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v)
45 {
46   PetscErrorCode   ierr;
47   PetscBool        flg = PETSC_FALSE, set;
48   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;
49 
50   PetscFunctionBegin;
51   ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr);
52   ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr);
53   ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr);
54   ierr = PetscOptionsBool("-viewer_hdf5_collective","Enable collective transfer mode","PetscViewerHDF5SetCollective",flg,&flg,&set);CHKERRQ(ierr);
55   if (set) {ierr = PetscViewerHDF5SetCollective(v,flg);CHKERRQ(ierr);}
56   ierr = PetscOptionsTail();CHKERRQ(ierr);
57   PetscFunctionReturn(0);
58 }
59 
60 static PetscErrorCode PetscViewerView_HDF5(PetscViewer v,PetscViewer viewer)
61 {
62   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*)v->data;
63   PetscBool         flg;
64   PetscErrorCode    ierr;
65 
66   PetscFunctionBegin;
67   if (hdf5->filename) {
68     ierr = PetscViewerASCIIPrintf(viewer,"Filename: %s\n",hdf5->filename);CHKERRQ(ierr);
69   }
70   ierr = PetscViewerASCIIPrintf(viewer,"Vectors with blocksize 1 saved as 2D datasets: %s\n",PetscBools[hdf5->basedimension2]);CHKERRQ(ierr);
71   ierr = PetscViewerASCIIPrintf(viewer,"Enforce single precision storage: %s\n",PetscBools[hdf5->spoutput]);CHKERRQ(ierr);
72   ierr = PetscViewerHDF5GetCollective(v,&flg);CHKERRQ(ierr);
73   ierr = PetscViewerASCIIPrintf(viewer,"MPI-IO transfer mode: %s\n",flg ? "collective" : "independent");CHKERRQ(ierr);
74   PetscFunctionReturn(0);
75 }
76 
77 static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
78 {
79   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
80   PetscErrorCode   ierr;
81 
82   PetscFunctionBegin;
83   ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);
84   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
85   PetscFunctionReturn(0);
86 }
87 
88 static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
89 {
90   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
91   PetscErrorCode   ierr;
92 
93   PetscFunctionBegin;
94   PetscStackCallHDF5(H5Pclose,(hdf5->dxpl_id));
95   ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr);
96   while (hdf5->groups) {
97     PetscViewerHDF5GroupList *tmp = hdf5->groups->next;
98 
99     ierr         = PetscFree(hdf5->groups->name);CHKERRQ(ierr);
100     ierr         = PetscFree(hdf5->groups);CHKERRQ(ierr);
101     hdf5->groups = tmp;
102   }
103   ierr = PetscFree(hdf5);CHKERRQ(ierr);
104   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
105   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr);
106   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
107   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr);
108   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr);
109   PetscFunctionReturn(0);
110 }
111 
112 static PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
113 {
114   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
115 
116   PetscFunctionBegin;
117   hdf5->btype = type;
118   PetscFunctionReturn(0);
119 }
120 
121 static PetscErrorCode  PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
122 {
123   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
124 
125   PetscFunctionBegin;
126   *type = hdf5->btype;
127   PetscFunctionReturn(0);
128 }
129 
130 static PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
131 {
132   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
133 
134   PetscFunctionBegin;
135   hdf5->basedimension2 = flg;
136   PetscFunctionReturn(0);
137 }
138 
139 /*@
140      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
141        dimension of 2.
142 
143     Logically Collective on PetscViewer
144 
145   Input Parameters:
146 +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
147 -  flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1
148 
149   Options Database:
150 .  -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
151 
152 
153   Notes:
154     Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
155          of one when the dimension is lower. Others think the option is crazy.
156 
157   Level: intermediate
158 
159 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
160 
161 @*/
162 PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
163 {
164   PetscErrorCode ierr;
165 
166   PetscFunctionBegin;
167   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
168   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
169   PetscFunctionReturn(0);
170 }
171 
172 /*@
173      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
174        dimension of 2.
175 
176     Logically Collective on PetscViewer
177 
178   Input Parameter:
179 .  viewer - the PetscViewer, must be of type HDF5
180 
181   Output Parameter:
182 .  flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1
183 
184   Notes:
185     Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
186          of one when the dimension is lower. Others think the option is crazy.
187 
188   Level: intermediate
189 
190 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
191 
192 @*/
193 PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
194 {
195   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
196 
197   PetscFunctionBegin;
198   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
199   *flg = hdf5->basedimension2;
200   PetscFunctionReturn(0);
201 }
202 
203 static PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
204 {
205   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
206 
207   PetscFunctionBegin;
208   hdf5->spoutput = flg;
209   PetscFunctionReturn(0);
210 }
211 
212 /*@
213      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
214        compiled with double precision PetscReal.
215 
216     Logically Collective on PetscViewer
217 
218   Input Parameters:
219 +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
220 -  flg - if PETSC_TRUE the data will be written to disk with single precision
221 
222   Options Database:
223 .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
224 
225 
226   Notes:
227     Setting this option does not make any difference if PETSc is compiled with single precision
228          in the first place. It does not affect reading datasets (HDF5 handle this internally).
229 
230   Level: intermediate
231 
232 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
233           PetscReal
234 
235 @*/
236 PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
237 {
238   PetscErrorCode ierr;
239 
240   PetscFunctionBegin;
241   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
242   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
243   PetscFunctionReturn(0);
244 }
245 
246 /*@
247      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
248        compiled with double precision PetscReal.
249 
250     Logically Collective on PetscViewer
251 
252   Input Parameter:
253 .  viewer - the PetscViewer, must be of type HDF5
254 
255   Output Parameter:
256 .  flg - if PETSC_TRUE the data will be written to disk with single precision
257 
258   Notes:
259     Setting this option does not make any difference if PETSc is compiled with single precision
260          in the first place. It does not affect reading datasets (HDF5 handle this internally).
261 
262   Level: intermediate
263 
264 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
265           PetscReal
266 
267 @*/
268 PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
269 {
270   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
271 
272   PetscFunctionBegin;
273   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
274   *flg = hdf5->spoutput;
275   PetscFunctionReturn(0);
276 }
277 
278 static PetscErrorCode  PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg)
279 {
280   PetscFunctionBegin;
281   /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions
282      - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */
283 #if H5_VERSION_GE(1,10,3) && defined(H5_HAVE_PARALLEL)
284   {
285     PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
286     PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT));
287   }
288 #else
289   if (flg) {
290     PetscErrorCode ierr;
291     ierr = PetscPrintf(PetscObjectComm((PetscObject)viewer), "Warning: PetscViewerHDF5SetCollective(viewer,PETSC_TRUE) is ignored for HDF5 versions prior to 1.10.3 or if built without MPI support\n");CHKERRQ(ierr);
292   }
293 #endif
294   PetscFunctionReturn(0);
295 }
296 
297 /*@
298   PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes.
299 
300   Logically Collective; flg must contain common value
301 
302   Input Parameters:
303 + viewer - the PetscViewer; if it is not hdf5 then this command is ignored
304 - flg - PETSC_TRUE for collective mode; PETSC_FALSE for independent mode (default)
305 
306   Options Database:
307 . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers
308 
309   Notes:
310   Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better.
311   However, this works correctly only since HDF5 1.10.3 and if HDF5 is installed for MPI; hence, we ignore this setting for older versions.
312 
313   Developer notes:
314   In the HDF5 layer, PETSC_TRUE / PETSC_FALSE means H5Pset_dxpl_mpio() is called with H5FD_MPIO_COLLECTIVE / H5FD_MPIO_INDEPENDENT, respectively.
315   This in turn means use of MPI_File_{read,write}_all /  MPI_File_{read,write} in the MPI-IO layer, respectively.
316   See HDF5 documentation and MPI-IO documentation for details.
317 
318   Level: intermediate
319 
320 .seealso: PetscViewerHDF5GetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open()
321 
322 @*/
323 PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg)
324 {
325   PetscErrorCode ierr;
326 
327   PetscFunctionBegin;
328   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
329   PetscValidLogicalCollectiveBool(viewer,flg,2);
330   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
331   PetscFunctionReturn(0);
332 }
333 
334 static PetscErrorCode  PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg)
335 {
336 #if defined(H5_HAVE_PARALLEL)
337   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*) viewer->data;
338   H5FD_mpio_xfer_t  mode;
339 #endif
340 
341   PetscFunctionBegin;
342 #if !defined(H5_HAVE_PARALLEL)
343   *flg = PETSC_FALSE;
344 #else
345   PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode));
346   *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE;
347 #endif
348   PetscFunctionReturn(0);
349 }
350 
351 /*@
352   PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes.
353 
354   Not Collective
355 
356   Input Parameters:
357 . viewer - the HDF5 PetscViewer
358 
359   Output Parameters:
360 . flg - the flag
361 
362   Level: intermediate
363 
364   Notes:
365   This setting works correctly only since HDF5 1.10.3 and if HDF5 was installed for MPI. For older versions, PETSC_FALSE will be always returned.
366   For more details, see PetscViewerHDF5SetCollective().
367 
368 .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open()
369 
370 @*/
371 PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg)
372 {
373   PetscErrorCode ierr;
374 
375   PetscFunctionBegin;
376   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
377   PetscValidBoolPointer(flg,2);
378 
379   ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));CHKERRQ(ierr);
380   PetscFunctionReturn(0);
381 }
382 
383 static PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
384 {
385   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
386   hid_t             plist_id;
387   PetscErrorCode    ierr;
388 
389   PetscFunctionBegin;
390   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
391   if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);}
392   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
393   /* Set up file access property list with parallel I/O access */
394   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
395 #if defined(H5_HAVE_PARALLEL)
396   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL));
397 #endif
398   /* Create or open the file collectively */
399   switch (hdf5->btype) {
400   case FILE_MODE_READ:
401     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
402     break;
403   case FILE_MODE_APPEND:
404   case FILE_MODE_UPDATE:
405   {
406     PetscBool flg;
407     ierr = PetscTestFile(hdf5->filename, 'r', &flg);CHKERRQ(ierr);
408     if (flg) PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
409     else     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id));
410     break;
411   }
412   case FILE_MODE_WRITE:
413     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
414     break;
415   case FILE_MODE_UNDEFINED:
416     SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
417   default:
418     SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP, "Unsupported file mode %s",PetscFileModes[hdf5->btype]);
419   }
420   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
421   PetscStackCallHDF5(H5Pclose,(plist_id));
422   PetscFunctionReturn(0);
423 }
424 
425 static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
426 {
427   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
428 
429   PetscFunctionBegin;
430   *name = vhdf5->filename;
431   PetscFunctionReturn(0);
432 }
433 
434 static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
435 {
436   /*
437   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
438   PetscErrorCode   ierr;
439   */
440 
441   PetscFunctionBegin;
442   PetscFunctionReturn(0);
443 }
444 
445 /*MC
446    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
447 
448 
449 .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
450            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
451            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
452            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
453 
454   Level: beginner
455 M*/
456 
457 PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
458 {
459   PetscViewer_HDF5 *hdf5;
460   PetscErrorCode   ierr;
461 
462   PetscFunctionBegin;
463   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
464 
465   v->data                = (void*) hdf5;
466   v->ops->destroy        = PetscViewerDestroy_HDF5;
467   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
468   v->ops->setup          = PetscViewerSetUp_HDF5;
469   v->ops->view           = PetscViewerView_HDF5;
470   v->ops->flush          = NULL;
471   hdf5->btype            = FILE_MODE_UNDEFINED;
472   hdf5->filename         = NULL;
473   hdf5->timestep         = -1;
474   hdf5->groups           = NULL;
475 
476   PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER));
477 
478   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
479   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr);
480   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
481   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr);
482   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
483   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr);
484   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);CHKERRQ(ierr);
485   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);CHKERRQ(ierr);
486   PetscFunctionReturn(0);
487 }
488 
489 /*@C
490    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
491 
492    Collective
493 
494    Input Parameters:
495 +  comm - MPI communicator
496 .  name - name of file
497 -  type - type of file
498 
499    Output Parameter:
500 .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
501 
502   Options Database:
503 +  -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
504 -  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
505 
506    Level: beginner
507 
508    Notes:
509    Reading is always available, regardless of the mode. Available modes are
510 +  FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY]
511 .  FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC]
512 .  FILE_MODE_APPEND - if file exists, keep existing contents [H5Fopen() with H5F_ACC_RDWR], else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_EXCL]
513 -  FILE_MODE_UPDATE - same as FILE_MODE_APPEND
514 
515    In case of FILE_MODE_APPEND / FILE_MODE_UPDATE, any stored object (dataset, attribute) can be selectively ovewritten if the same fully qualified name (/group/path/to/object) is specified.
516 
517    This PetscViewer should be destroyed with PetscViewerDestroy().
518 
519 
520 .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
521           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
522           MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
523 @*/
524 PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
525 {
526   PetscErrorCode ierr;
527 
528   PetscFunctionBegin;
529   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
530   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
531   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
532   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
533   ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr);
534   PetscFunctionReturn(0);
535 }
536 
537 /*@C
538   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
539 
540   Not collective
541 
542   Input Parameter:
543 . viewer - the PetscViewer
544 
545   Output Parameter:
546 . file_id - The file id
547 
548   Level: intermediate
549 
550 .seealso: PetscViewerHDF5Open()
551 @*/
552 PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
553 {
554   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
555 
556   PetscFunctionBegin;
557   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
558   if (file_id) *file_id = hdf5->file_id;
559   PetscFunctionReturn(0);
560 }
561 
562 /*@C
563   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
564 
565   Not collective
566 
567   Input Parameters:
568 + viewer - the PetscViewer
569 - name - The group name
570 
571   Level: intermediate
572 
573   Notes:
574   This is designed to mnemonically resemble the Unix cd command.
575   + If name begins with '/', it is interpreted as an absolute path fully replacing current group, otherwise it is taken as relative to the current group.
576   . NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/".
577   - "." means the current group is pushed again.
578 
579   Example:
580   Suppose the current group is "/a".
581   + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/".
582   . If name is ".", then the new group will be "/a".
583   . If name is "b", then the new group will be "/a/b".
584   - If name is "/b", then the new group will be "/b".
585 
586   Developer Notes:
587   The root group "/" is internally stored as NULL.
588 
589 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
590 @*/
591 PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[])
592 {
593   PetscViewer_HDF5          *hdf5 = (PetscViewer_HDF5*) viewer->data;
594   PetscViewerHDF5GroupList  *groupNode;
595   size_t                    i,len;
596   char                      buf[PETSC_MAX_PATH_LEN];
597   const char                *gname;
598   PetscErrorCode            ierr;
599 
600   PetscFunctionBegin;
601   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
602   if (name) PetscValidCharPointer(name,2);
603   ierr = PetscStrlen(name, &len);CHKERRQ(ierr);
604   gname = NULL;
605   if (len) {
606     if (len == 1 && name[0] == '.') {
607       /* use current name */
608       gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL;
609     } else if (name[0] == '/') {
610       /* absolute */
611       for (i=1; i<len; i++) {
612         if (name[i] != '/') {
613           gname = name;
614           break;
615         }
616       }
617     } else {
618       /* relative */
619       const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : "";
620       ierr = PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name);CHKERRQ(ierr);
621       gname = buf;
622     }
623   }
624   ierr = PetscNew(&groupNode);CHKERRQ(ierr);
625   ierr = PetscStrallocpy(gname, (char**) &groupNode->name);CHKERRQ(ierr);
626   groupNode->next = hdf5->groups;
627   hdf5->groups    = groupNode;
628   PetscFunctionReturn(0);
629 }
630 
631 /*@
632   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
633 
634   Not collective
635 
636   Input Parameter:
637 . viewer - the PetscViewer
638 
639   Level: intermediate
640 
641 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
642 @*/
643 PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
644 {
645   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
646   PetscViewerHDF5GroupList *groupNode;
647   PetscErrorCode   ierr;
648 
649   PetscFunctionBegin;
650   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
651   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
652   groupNode    = hdf5->groups;
653   hdf5->groups = hdf5->groups->next;
654   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
655   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
656   PetscFunctionReturn(0);
657 }
658 
659 /*@C
660   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup().
661   If none has been assigned, returns NULL.
662 
663   Not collective
664 
665   Input Parameter:
666 . viewer - the PetscViewer
667 
668   Output Parameter:
669 . name - The group name
670 
671   Level: intermediate
672 
673 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
674 @*/
675 PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[])
676 {
677   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
678 
679   PetscFunctionBegin;
680   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
681   PetscValidPointer(name,2);
682   if (hdf5->groups) *name = hdf5->groups->name;
683   else *name = NULL;
684   PetscFunctionReturn(0);
685 }
686 
687 /*@
688   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(),
689   and return this group's ID and file ID.
690   If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID.
691 
692   Not collective
693 
694   Input Parameter:
695 . viewer - the PetscViewer
696 
697   Output Parameter:
698 + fileId - The HDF5 file ID
699 - groupId - The HDF5 group ID
700 
701   Notes:
702   If the viewer is writable, the group is created if it doesn't exist yet.
703 
704   Level: intermediate
705 
706 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
707 @*/
708 PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
709 {
710   hid_t          file_id;
711   H5O_type_t     type;
712   const char     *groupName = NULL;
713   PetscBool      create;
714   PetscErrorCode ierr;
715 
716   PetscFunctionBegin;
717   ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr);
718   ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr);
719   ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr);
720   ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr);
721   if (type != H5O_TYPE_GROUP) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s resolves to something which is not a group", groupName);
722   PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT));
723   *fileId  = file_id;
724   PetscFunctionReturn(0);
725 }
726 
727 /*@
728   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
729 
730   Not collective
731 
732   Input Parameter:
733 . viewer - the PetscViewer
734 
735   Level: intermediate
736 
737 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
738 @*/
739 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
740 {
741   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
742 
743   PetscFunctionBegin;
744   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
745   ++hdf5->timestep;
746   PetscFunctionReturn(0);
747 }
748 
749 /*@
750   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
751   of -1 disables blocking with timesteps.
752 
753   Not collective
754 
755   Input Parameters:
756 + viewer - the PetscViewer
757 - timestep - The timestep number
758 
759   Level: intermediate
760 
761 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
762 @*/
763 PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
764 {
765   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
766 
767   PetscFunctionBegin;
768   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
769   hdf5->timestep = timestep;
770   PetscFunctionReturn(0);
771 }
772 
773 /*@
774   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
775 
776   Not collective
777 
778   Input Parameter:
779 . viewer - the PetscViewer
780 
781   Output Parameter:
782 . timestep - The timestep number
783 
784   Level: intermediate
785 
786 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
787 @*/
788 PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
789 {
790   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
791 
792   PetscFunctionBegin;
793   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
794   PetscValidPointer(timestep,2);
795   *timestep = hdf5->timestep;
796   PetscFunctionReturn(0);
797 }
798 
799 /*@C
800   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
801 
802   Not collective
803 
804   Input Parameter:
805 . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
806 
807   Output Parameter:
808 . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
809 
810   Level: advanced
811 
812 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
813 @*/
814 PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
815 {
816   PetscFunctionBegin;
817   if (ptype == PETSC_INT)
818 #if defined(PETSC_USE_64BIT_INDICES)
819                                        *htype = H5T_NATIVE_LLONG;
820 #else
821                                        *htype = H5T_NATIVE_INT;
822 #endif
823   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
824   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
825   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
826   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_INT;
827   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_INT;
828   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
829   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
830   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
831   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
832   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
833   PetscFunctionReturn(0);
834 }
835 
836 /*@C
837   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
838 
839   Not collective
840 
841   Input Parameter:
842 . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
843 
844   Output Parameter:
845 . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
846 
847   Level: advanced
848 
849 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
850 @*/
851 PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
852 {
853   PetscFunctionBegin;
854 #if defined(PETSC_USE_64BIT_INDICES)
855   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
856   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
857 #else
858   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
859 #endif
860   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
861   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
862   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
863   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
864   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
865   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
866   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
867   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
868   PetscFunctionReturn(0);
869 }
870 
871 /*@C
872  PetscViewerHDF5WriteAttribute - Write an attribute
873 
874   Input Parameters:
875 + viewer - The HDF5 viewer
876 . parent - The parent dataset/group name
877 . name   - The attribute name
878 . datatype - The attribute type
879 - value    - The attribute value
880 
881   Level: advanced
882 
883   Notes:
884   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else the dataset/group is relative to the current pushed group. NULL means the current pushed group.
885 
886 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
887 @*/
888 PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
889 {
890   char           *parentAbsPath;
891   hid_t          h5, dataspace, obj, attribute, dtype;
892   PetscBool      has;
893   PetscErrorCode ierr;
894 
895   PetscFunctionBegin;
896   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
897   if (parent) PetscValidCharPointer(parent, 2);
898   PetscValidCharPointer(name, 3);
899   PetscValidLogicalCollectiveEnum(viewer,datatype,4);
900   PetscValidPointer(value, 5);
901   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr);
902   ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr);
903   ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);CHKERRQ(ierr);
904   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
905   if (datatype == PETSC_STRING) {
906     size_t len;
907     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
908     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
909   }
910   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
911   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
912   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT));
913   if (has) {
914     PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
915   } else {
916     PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
917   }
918   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
919   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
920   PetscStackCallHDF5(H5Aclose,(attribute));
921   PetscStackCallHDF5(H5Oclose,(obj));
922   PetscStackCallHDF5(H5Sclose,(dataspace));
923   ierr = PetscFree(parentAbsPath);CHKERRQ(ierr);
924   PetscFunctionReturn(0);
925 }
926 
927 /*@C
928  PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name
929 
930   Input Parameters:
931 + viewer   - The HDF5 viewer
932 . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
933 . name     - The attribute name
934 . datatype - The attribute type
935 - value    - The attribute value
936 
937   Notes:
938   This fails if currentGroup/objectName doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
939   You might want to check first if it does using PetscViewerHDF5HasObject().
940 
941   Level: advanced
942 
943 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
944 @*/
945 PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
946 {
947   PetscErrorCode ierr;
948 
949   PetscFunctionBegin;
950   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
951   PetscValidHeader(obj,2);
952   PetscValidCharPointer(name,3);
953   PetscValidPointer(value,5);
954   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
955   ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
956   PetscFunctionReturn(0);
957 }
958 
959 /*@C
960  PetscViewerHDF5ReadAttribute - Read an attribute
961 
962   Input Parameters:
963 + viewer - The HDF5 viewer
964 . parent - The parent dataset/group name
965 . name   - The attribute name
966 - datatype - The attribute type
967 
968   Output Parameter:
969 . value    - The attribute value
970 
971   Notes: If the datatype is PETSC_STRING one must PetscFree() the obtained value when it is no longer needed.
972 
973   Level: advanced
974 
975   Notes:
976   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else the dataset/group is relative to the current pushed group. NULL means the current pushed group.
977 
978 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
979 @*/
980 PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, void *value)
981 {
982   char           *parentAbsPath;
983   hid_t          h5, obj, attribute, atype, dtype;
984   PetscBool      has;
985   PetscErrorCode ierr;
986 
987   PetscFunctionBegin;
988   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
989   if (parent) PetscValidCharPointer(parent, 2);
990   PetscValidCharPointer(name, 3);
991   PetscValidPointer(value, 5);
992   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr);
993   ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL);CHKERRQ(ierr);
994   if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);CHKERRQ(ierr);}
995   if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parentAbsPath, name);
996   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
997   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
998   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT));
999   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
1000   if (datatype == PETSC_STRING) {
1001     size_t len;
1002     PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
1003     PetscStackCall("H5Tget_size",len = H5Tget_size(atype));
1004     ierr = PetscMalloc((len+1) * sizeof(char), value);CHKERRQ(ierr);
1005     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
1006     PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value));
1007   } else {
1008     PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
1009   }
1010   PetscStackCallHDF5(H5Aclose,(attribute));
1011   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
1012   PetscStackCallHDF5(H5Oclose,(obj));
1013   ierr = PetscFree(parentAbsPath);CHKERRQ(ierr);
1014   PetscFunctionReturn(0);
1015 }
1016 
1017 /*@C
1018  PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name
1019 
1020   Input Parameters:
1021 + viewer   - The HDF5 viewer
1022 . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1023 . name     - The attribute name
1024 - datatype - The attribute type
1025 
1026   Output Parameter:
1027 . value    - The attribute value
1028 
1029   Notes:
1030   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1031   You might want to check first if it does using PetscViewerHDF5HasObject().
1032 
1033   Level: advanced
1034 
1035 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1036 @*/
1037 PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value)
1038 {
1039   PetscErrorCode ierr;
1040 
1041   PetscFunctionBegin;
1042   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1043   PetscValidHeader(obj,2);
1044   PetscValidCharPointer(name,3);
1045   PetscValidPointer(value, 5);
1046   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
1047   ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
1048   PetscFunctionReturn(0);
1049 }
1050 
1051 PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
1052 {
1053   htri_t exists;
1054   hid_t group;
1055 
1056   PetscFunctionBegin;
1057   PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT));
1058   if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT));
1059   if (!exists && createGroup) {
1060     PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1061     PetscStackCallHDF5(H5Gclose,(group));
1062     exists = PETSC_TRUE;
1063   }
1064   *exists_ = (PetscBool) exists;
1065   PetscFunctionReturn(0);
1066 }
1067 
1068 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
1069 {
1070   const char     rootGroupName[] = "/";
1071   hid_t          h5;
1072   PetscBool      exists=PETSC_FALSE;
1073   PetscInt       i;
1074   int            n;
1075   char           **hierarchy;
1076   char           buf[PETSC_MAX_PATH_LEN]="";
1077   PetscErrorCode ierr;
1078 
1079   PetscFunctionBegin;
1080   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1081   if (name) PetscValidCharPointer(name, 2);
1082   else name = rootGroupName;
1083   if (has) {
1084     PetscValidIntPointer(has, 3);
1085     *has = PETSC_FALSE;
1086   }
1087   if (otype) {
1088     PetscValidIntPointer(otype, 4);
1089     *otype = H5O_TYPE_UNKNOWN;
1090   }
1091   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
1092 
1093   /*
1094      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
1095      Hence, each of them needs to be tested separately:
1096      1) whether it's a valid link
1097      2) whether this link resolves to an object
1098      See H5Oexists_by_name() documentation.
1099   */
1100   ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr);
1101   if (!n) {
1102     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1103     if (has)   *has   = PETSC_TRUE;
1104     if (otype) *otype = H5O_TYPE_GROUP;
1105     ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
1106     PetscFunctionReturn(0);
1107   }
1108   for (i=0; i<n; i++) {
1109     ierr = PetscStrcat(buf,"/");CHKERRQ(ierr);
1110     ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr);
1111     ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr);
1112     if (!exists) break;
1113   }
1114   ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
1115 
1116   /* If the object exists, get its type */
1117   if (exists && otype) {
1118     H5O_info_t info;
1119 
1120     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
1121     PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT));
1122     *otype = info.type;
1123   }
1124   if (has) *has = exists;
1125   PetscFunctionReturn(0);
1126 }
1127 
1128 /*@C
1129  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
1130 
1131   Input Parameters:
1132 + viewer - The HDF5 viewer
1133 - path - The group path
1134 
1135   Output Parameter:
1136 . has - Flag for group existence
1137 
1138   Level: advanced
1139 
1140   Notes:
1141   If path starts with '/', it is taken as an absolute path overriding currently pushed group, else the dataset/group is relative to the current pushed group. NULL means the current pushed group.
1142   If the path exists but is not a group, PETSC_FALSE is returned.
1143 
1144 .seealso: PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasDataset(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup(), PetscViewerHDF5OpenGroup()
1145 @*/
1146 PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has)
1147 {
1148   H5O_type_t     type;
1149   char           *abspath;
1150   PetscErrorCode ierr;
1151 
1152   PetscFunctionBegin;
1153   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1154   if (path) PetscValidCharPointer(path,2);
1155   PetscValidBoolPointer(has,3);
1156   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);CHKERRQ(ierr);
1157   ierr = PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);CHKERRQ(ierr);
1158   *has = (PetscBool)(type == H5O_TYPE_GROUP);
1159   ierr = PetscFree(abspath);CHKERRQ(ierr);
1160   PetscFunctionReturn(0);
1161 }
1162 
1163 /*@C
1164  PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file
1165 
1166   Input Parameters:
1167 + viewer - The HDF5 viewer
1168 - path - The dataset path
1169 
1170   Output Parameter:
1171 . has - Flag whether dataset exists
1172 
1173   Level: advanced
1174 
1175   Notes:
1176   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else the dataset/group is relative to the current pushed group.
1177   If path is NULL or empty, has is set to PETSC_FALSE.
1178   If the path exists but is not a dataset, has is set to PETSC_FALSE as well.
1179 
1180 .seealso: PetscViewerHDF5HasObject(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasGroup(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup()
1181 @*/
1182 PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has)
1183 {
1184   H5O_type_t     type;
1185   char           *abspath;
1186   PetscErrorCode ierr;
1187 
1188   PetscFunctionBegin;
1189   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1190   if (path) PetscValidCharPointer(path,2);
1191   PetscValidBoolPointer(has,3);
1192   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);CHKERRQ(ierr);
1193   ierr = PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);CHKERRQ(ierr);
1194   *has = (PetscBool)(type == H5O_TYPE_DATASET);
1195   ierr = PetscFree(abspath);CHKERRQ(ierr);
1196   PetscFunctionReturn(0);
1197 }
1198 
1199 /*@
1200  PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1201 
1202   Input Parameters:
1203 + viewer - The HDF5 viewer
1204 - obj    - The named object
1205 
1206   Output Parameter:
1207 . has    - Flag for dataset existence
1208 
1209   Notes:
1210   If the object is unnamed, an error occurs.
1211   If the path currentGroup/objectName exists but is not a dataset, has is set to PETSC_FALSE as well.
1212 
1213   Level: advanced
1214 
1215 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasDataset(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup()
1216 @*/
1217 PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1218 {
1219   size_t         len;
1220   PetscErrorCode ierr;
1221 
1222   PetscFunctionBegin;
1223   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1224   PetscValidHeader(obj,2);
1225   PetscValidBoolPointer(has,3);
1226   ierr = PetscStrlen(obj->name, &len);CHKERRQ(ierr);
1227   if (!len) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
1228   ierr = PetscViewerHDF5HasDataset(viewer, obj->name, has);CHKERRQ(ierr);
1229   PetscFunctionReturn(0);
1230 }
1231 
1232 /*@C
1233  PetscViewerHDF5HasAttribute - Check whether an attribute exists
1234 
1235   Input Parameters:
1236 + viewer - The HDF5 viewer
1237 . parent - The parent dataset/group name
1238 - name   - The attribute name
1239 
1240   Output Parameter:
1241 . has    - Flag for attribute existence
1242 
1243   Level: advanced
1244 
1245   Notes:
1246   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else the dataset/group is relative to the current pushed group. NULL means the current pushed group.
1247 
1248 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1249 @*/
1250 PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1251 {
1252   char           *parentAbsPath;
1253   PetscErrorCode ierr;
1254 
1255   PetscFunctionBegin;
1256   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1257   if (parent) PetscValidCharPointer(parent,2);
1258   PetscValidCharPointer(name,3);
1259   PetscValidBoolPointer(has,4);
1260   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr);
1261   ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL);CHKERRQ(ierr);
1262   if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has);CHKERRQ(ierr);}
1263   ierr = PetscFree(parentAbsPath);CHKERRQ(ierr);
1264   PetscFunctionReturn(0);
1265 }
1266 
1267 /*@C
1268  PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name
1269 
1270   Input Parameters:
1271 + viewer - The HDF5 viewer
1272 . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1273 - name   - The attribute name
1274 
1275   Output Parameter:
1276 . has    - Flag for attribute existence
1277 
1278   Notes:
1279   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1280   You might want to check first if it does using PetscViewerHDF5HasObject().
1281 
1282   Level: advanced
1283 
1284 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1285 @*/
1286 PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1287 {
1288   PetscErrorCode ierr;
1289 
1290   PetscFunctionBegin;
1291   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1292   PetscValidHeader(obj,2);
1293   PetscValidCharPointer(name,3);
1294   PetscValidBoolPointer(has,4);
1295   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
1296   ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr);
1297   PetscFunctionReturn(0);
1298 }
1299 
1300 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1301 {
1302   hid_t          h5;
1303   htri_t         hhas;
1304   PetscErrorCode ierr;
1305 
1306   PetscFunctionBegin;
1307   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
1308   PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT));
1309   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1310   PetscFunctionReturn(0);
1311 }
1312 
1313 /*
1314   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1315   is attached to a communicator, in this case the attribute is a PetscViewer.
1316 */
1317 PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1318 
1319 /*@C
1320   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
1321 
1322   Collective
1323 
1324   Input Parameter:
1325 . comm - the MPI communicator to share the HDF5 PetscViewer
1326 
1327   Level: intermediate
1328 
1329   Options Database Keys:
1330 . -viewer_hdf5_filename <name>
1331 
1332   Environmental variables:
1333 . PETSC_VIEWER_HDF5_FILENAME
1334 
1335   Notes:
1336   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
1337   an error code.  The HDF5 PetscViewer is usually used in the form
1338 $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1339 
1340 .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1341 @*/
1342 PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1343 {
1344   PetscErrorCode ierr;
1345   PetscBool      flg;
1346   PetscViewer    viewer;
1347   char           fname[PETSC_MAX_PATH_LEN];
1348   MPI_Comm       ncomm;
1349 
1350   PetscFunctionBegin;
1351   ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1352   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1353     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,NULL);
1354     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1355   }
1356   ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1357   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1358   if (!flg) { /* PetscViewer not yet created */
1359     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
1360     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1361     if (!flg) {
1362       ierr = PetscStrcpy(fname,"output.h5");
1363       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1364     }
1365     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
1366     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1367     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
1368     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1369     ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
1370     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1371   }
1372   ierr = PetscCommDestroy(&ncomm);
1373   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1374   PetscFunctionReturn(viewer);
1375 }
1376