Actual source code: hdf5v.c

  1: #include <petsc/private/petscimpl.h>
  2: #include <petscviewerhdf5.h>

  4: /*@C
  5:   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with `PetscViewerHDF5PushGroup()`/`PetscViewerHDF5PopGroup()`.

  7:   Not Collective

  9:   Input Parameters:
 10: + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
 11: - path   - (Optional) The path relative to the pushed group

 13:   Output Parameter:
 14: . abspath - The absolute HDF5 path (group)

 16:   Level: intermediate

 18:   Notes:
 19:   If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group.
 20:   `NULL` or empty path means the current pushed group.

 22:   The output abspath is newly allocated so needs to be freed.

 24: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()`
 25: @*/
 26: PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char path[], const char *abspath[])
 27: {
 28:   PetscFunctionBegin;
 30:   if (path) PetscAssertPointer(path, 2);
 31:   PetscAssertPointer(abspath, 3);
 32:   PetscUseMethod(viewer, "PetscViewerHDF5GetGroup_C", (PetscViewer, const char[], const char *[]), (viewer, path, abspath));
 33:   PetscFunctionReturn(PETSC_SUCCESS);
 34: }

 36: static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj)
 37: {
 38:   PetscBool has;

 40:   PetscFunctionBegin;
 41:   PetscCall(PetscViewerHDF5HasObject(viewer, obj, &has));
 42:   if (!has) {
 43:     const char *group;
 44:     PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group));
 45:     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", obj->name, group);
 46:   }
 47:   PetscFunctionReturn(PETSC_SUCCESS);
 48: }

 50: /*@
 51:   PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
 52:   dimension of 2.

 54:   Logically Collective

 56:   Input Parameters:
 57: + viewer - the `PetscViewer`; if it is a `PETSCVIEWERHDF5` then this command is ignored
 58: - flg    - if `PETSC_TRUE` the vector will always have at least a dimension of 2 even if that first dimension is of size 1

 60:   Options Database Key:
 61: . -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

 63:   Level: intermediate

 65:   Note:
 66:   Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
 67:   of one when the dimension is lower. Others think the option is crazy.

 69: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
 70: @*/
 71: PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer, PetscBool flg)
 72: {
 73:   PetscFunctionBegin;
 75:   PetscTryMethod(viewer, "PetscViewerHDF5SetBaseDimension2_C", (PetscViewer, PetscBool), (viewer, flg));
 76:   PetscFunctionReturn(PETSC_SUCCESS);
 77: }

 79: /*@
 80:   PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
 81:   dimension of 2.

 83:   Logically Collective

 85:   Input Parameter:
 86: . viewer - the `PetscViewer`, must be `PETSCVIEWERHDF5`

 88:   Output Parameter:
 89: . flg - if `PETSC_TRUE` the vector will always have at least a dimension of 2 even if that first dimension is of size 1

 91:   Level: intermediate

 93:   Note:
 94:   Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
 95:   of one when the dimension is lower. Others think the option is crazy.

 97: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
 98: @*/
 99: PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer, PetscBool *flg)
100: {
101:   PetscFunctionBegin;
103:   PetscUseMethod(viewer, "PetscViewerHDF5GetBaseDimension2_C", (PetscViewer, PetscBool *), (viewer, flg));
104:   PetscFunctionReturn(PETSC_SUCCESS);
105: }

107: /*@
108:   PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
109:   compiled with double precision `PetscReal`.

111:   Logically Collective

113:   Input Parameters:
114: + viewer - the PetscViewer; if it is a `PETSCVIEWERHDF5` then this command is ignored
115: - flg    - if `PETSC_TRUE` the data will be written to disk with single precision

117:   Options Database Key:
118: . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision

120:   Level: intermediate

122:   Note:
123:   Setting this option does not make any difference if PETSc is compiled with single precision
124:   in the first place. It does not affect reading datasets (HDF5 handle this internally).

126: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
127:           `PetscReal`, `PetscViewerHDF5GetSPOutput()`
128: @*/
129: PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer, PetscBool flg)
130: {
131:   PetscFunctionBegin;
133:   PetscTryMethod(viewer, "PetscViewerHDF5SetSPOutput_C", (PetscViewer, PetscBool), (viewer, flg));
134:   PetscFunctionReturn(PETSC_SUCCESS);
135: }

137: /*@
138:   PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
139:   compiled with double precision `PetscReal`.

141:   Logically Collective

143:   Input Parameter:
144: . viewer - the PetscViewer, must be of type `PETSCVIEWERHDF5`

146:   Output Parameter:
147: . flg - if `PETSC_TRUE` the data will be written to disk with single precision

149:   Level: intermediate

151: .seealso: [](sec_viewers), `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
152:           `PetscReal`, `PetscViewerHDF5SetSPOutput()`
153: @*/
154: PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer, PetscBool *flg)
155: {
156:   PetscFunctionBegin;
158:   PetscUseMethod(viewer, "PetscViewerHDF5GetSPOutput_C", (PetscViewer, PetscBool *), (viewer, flg));
159:   PetscFunctionReturn(PETSC_SUCCESS);
160: }

162: /*@
163:   PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes.

165:   Logically Collective; flg must contain common value

167:   Input Parameters:
168: + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored
169: - flg    - `PETSC_TRUE` for collective mode; `PETSC_FALSE` for independent mode (default)

171:   Options Database Key:
172: . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers

174:   Level: intermediate

176:   Note:
177:   Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better.
178:   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.

180:   Developer Notes:
181:   In the HDF5 layer, `PETSC_TRUE` / `PETSC_FALSE` means `H5Pset_dxpl_mpio()` is called with `H5FD_MPIO_COLLECTIVE` / `H5FD_MPIO_INDEPENDENT`, respectively.
182:   This in turn means use of MPI_File_{read,write}_all /  MPI_File_{read,write} in the MPI-IO layer, respectively.
183:   See HDF5 documentation and MPI-IO documentation for details.

185: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5GetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()`
186: @*/
187: PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer, PetscBool flg)
188: {
189:   PetscFunctionBegin;
192:   PetscTryMethod(viewer, "PetscViewerHDF5SetCollective_C", (PetscViewer, PetscBool), (viewer, flg));
193:   PetscFunctionReturn(PETSC_SUCCESS);
194: }

196: /*@
197:   PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes.

199:   Not Collective

201:   Input Parameter:
202: . viewer - the `PETSCVIEWERHDF5` `PetscViewer`

204:   Output Parameter:
205: . flg - the flag

207:   Level: intermediate

209:   Note:
210:   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.
211:   For more details, see `PetscViewerHDF5SetCollective()`.

213: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5SetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()`
214: @*/
215: PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer, PetscBool *flg)
216: {
217:   PetscFunctionBegin;
219:   PetscAssertPointer(flg, 2);
220:   PetscUseMethod(viewer, "PetscViewerHDF5GetCollective_C", (PetscViewer, PetscBool *), (viewer, flg));
221:   PetscFunctionReturn(PETSC_SUCCESS);
222: }

224: /*@
225:   PetscViewerHDF5SetDefaultTimestepping - Set the flag for default timestepping

227:   Logically Collective

229:   Input Parameters:
230: + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored
231: - flg    - if `PETSC_TRUE` we will assume that timestepping is on

233:   Options Database Key:
234: . -viewer_hdf5_default_timestepping - turns on (true) or off (false) default timestepping

236:   Level: intermediate

238:   Note:
239:   If the timestepping attribute is not found for an object, then the default timestepping is used

241: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5GetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()`
242: @*/
243: PetscErrorCode PetscViewerHDF5SetDefaultTimestepping(PetscViewer viewer, PetscBool flg)
244: {
245:   PetscFunctionBegin;
247:   PetscTryMethod(viewer, "PetscViewerHDF5SetDefaultTimestepping_C", (PetscViewer, PetscBool), (viewer, flg));
248:   PetscFunctionReturn(PETSC_SUCCESS);
249: }

251: /*@
252:   PetscViewerHDF5GetDefaultTimestepping - Get the flag for default timestepping

254:   Not Collective

256:   Input Parameter:
257: . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`

259:   Output Parameter:
260: . flg - if `PETSC_TRUE` we will assume that timestepping is on

262:   Level: intermediate

264: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5SetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()`
265: @*/
266: PetscErrorCode PetscViewerHDF5GetDefaultTimestepping(PetscViewer viewer, PetscBool *flg)
267: {
268:   PetscFunctionBegin;
270:   PetscUseMethod(viewer, "PetscViewerHDF5GetDefaultTimestepping_C", (PetscViewer, PetscBool *), (viewer, flg));
271:   PetscFunctionReturn(PETSC_SUCCESS);
272: }

274: /*@
275:   PetscViewerHDF5SetCompress - Set the flag for compression

277:   Logically Collective

279:   Input Parameters:
280: + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored
281: - flg    - if `PETSC_TRUE` we will turn on compression

283:   Options Database Key:
284: . -viewer_hdf5_compress - turns on (true) or off (false) compression

286:   Level: intermediate

288: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5GetCompress()`
289: @*/
290: PetscErrorCode PetscViewerHDF5SetCompress(PetscViewer viewer, PetscBool flg)
291: {
292:   PetscFunctionBegin;
294:   PetscTryMethod(viewer, "PetscViewerHDF5SetCompress_C", (PetscViewer, PetscBool), (viewer, flg));
295:   PetscFunctionReturn(PETSC_SUCCESS);
296: }

298: /*@
299:   PetscViewerHDF5GetCompress - Get the flag for compression

301:   Not Collective

303:   Input Parameter:
304: . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`

306:   Output Parameter:
307: . flg - if `PETSC_TRUE` we will turn on compression

309:   Level: intermediate

311: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5SetCompress()`
312: @*/
313: PetscErrorCode PetscViewerHDF5GetCompress(PetscViewer viewer, PetscBool *flg)
314: {
315:   PetscFunctionBegin;
317:   PetscUseMethod(viewer, "PetscViewerHDF5GetCompress_C", (PetscViewer, PetscBool *), (viewer, flg));
318:   PetscFunctionReturn(PETSC_SUCCESS);
319: }

321: /*@C
322:   PetscViewerHDF5Open - Opens a file for HDF5 input/output as a `PETSCVIEWERHDF5` `PetscViewer`

324:   Collective

326:   Input Parameters:
327: + comm - MPI communicator
328: . name - name of file
329: - type - type of file

331:   Output Parameter:
332: . hdf5v - `PetscViewer` for HDF5 input/output to use with the specified file

334:   Options Database Keys:
335: + -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
336: - -viewer_hdf5_sp_output       - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal

338:   Level: beginner

340:   Notes:
341:   Reading is always available, regardless of the mode. Available modes are
342: .vb
343:   FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY]
344:   FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC]
345:   FILE_MODE_APPEND - if file exists, keep existing contents [H5Fopen() with H5F_ACC_RDWR], else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_EXCL]
346:   FILE_MODE_UPDATE - same as FILE_MODE_APPEND
347: .ve

349:   In case of `FILE_MODE_APPEND` / `FILE_MODE_UPDATE`, any stored object (dataset, attribute) can be selectively overwritten if the same fully qualified name (/group/path/to/object) is specified.

351:   This PetscViewer should be destroyed with PetscViewerDestroy().

353: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, `PetscViewerHDF5SetBaseDimension2()`,
354:           `PetscViewerHDF5SetSPOutput()`, `PetscViewerHDF5GetBaseDimension2()`, `VecView()`, `MatView()`, `VecLoad()`,
355:           `MatLoad()`, `PetscFileMode`, `PetscViewer`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()`
356: @*/
357: PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
358: {
359:   PetscFunctionBegin;
360:   PetscCall(PetscViewerCreate(comm, hdf5v));
361:   PetscCall(PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5));
362:   PetscCall(PetscViewerFileSetMode(*hdf5v, type));
363:   PetscCall(PetscViewerFileSetName(*hdf5v, name));
364:   PetscCall(PetscViewerSetFromOptions(*hdf5v));
365:   PetscFunctionReturn(PETSC_SUCCESS);
366: }

368: /*@C
369:   PetscViewerHDF5PushGroup - Set the current HDF5 group for output

371:   Not Collective

373:   Input Parameters:
374: + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
375: - name   - The group name

377:   Level: intermediate

379:   Notes:
380:   This is designed to mnemonically resemble the Unix cd command.
381: .vb
382:   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.
383:   `NULL`, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/".
384:   "." means the current group is pushed again.
385: .ve

387:   Example:
388:   Suppose the current group is "/a".
389: .vb
390:   If name is `NULL`, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/".
391:   If name is ".", then the new group will be "/a".
392:   If name is "b", then the new group will be "/a/b".
393:   If name is "/b", then the new group will be "/b".
394: .ve

396:   Developer Notes:
397:   The root group "/" is internally stored as `NULL`.

399: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()`
400: @*/
401: PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[])
402: {
403:   PetscFunctionBegin;
405:   if (name) PetscAssertPointer(name, 2);
406:   PetscUseMethod(viewer, "PetscViewerHDF5PushGroup_C", (PetscViewer, const char[]), (viewer, name));
407:   PetscFunctionReturn(PETSC_SUCCESS);
408: }

410: /*@
411:   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value

413:   Not Collective

415:   Input Parameter:
416: . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`

418:   Level: intermediate

420: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()`
421: @*/
422: PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer)
423: {
424:   PetscFunctionBegin;
426:   PetscUseMethod(viewer, "PetscViewerHDF5PopGroup_C", (PetscViewer), (viewer));
427:   PetscFunctionReturn(PETSC_SUCCESS);
428: }

430: /*@C
431:   PetscViewerHDF5WriteGroup - Ensure the HDF5 group exists in the HDF5 file

433:   Not Collective

435:   Input Parameters:
436: + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
437: - path   - (Optional) The path relative to the pushed group

439:   Level: intermediate

441:   Note:
442:   If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group.
443:   `NULL` or empty path means the current pushed group.

445:   This will fail if the viewer is not writable.

447: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
448: @*/
449: PetscErrorCode PetscViewerHDF5WriteGroup(PetscViewer viewer, const char path[])
450: {
451:   PetscFunctionBegin;
453:   PetscUseMethod(viewer, "PetscViewerHDF5WriteGroup_C", (PetscViewer, const char[]), (viewer, path));
454:   PetscFunctionReturn(PETSC_SUCCESS);
455: }

457: /*@
458:   PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing.

460:   Not Collective

462:   Input Parameter:
463: . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`

465:   Level: intermediate

467:   Notes:
468:   On first `PetscViewerHDF5PushTimestepping()`, the initial time step is set to 0.
469:   Next timesteps can then be set using `PetscViewerHDF5IncrementTimestep()` or `PetscViewerHDF5SetTimestep()`.
470:   Current timestep value determines which timestep is read from or written to any dataset on the next HDF5 I/O operation [e.g. `VecView()`].
471:   Use `PetscViewerHDF5PopTimestepping()` to deactivate timestepping mode; calling it by the end of the program is NOT mandatory.
472:   Current timestep is remembered between `PetscViewerHDF5PopTimestepping()` and the next `PetscViewerHDF5PushTimestepping()`.

474:   If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again.
475:   Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error.

477:   Developer Notes:
478:   Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true.

480: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
481: @*/
482: PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer)
483: {
484:   PetscFunctionBegin;
486:   PetscUseMethod(viewer, "PetscViewerHDF5PushTimestepping_C", (PetscViewer), (viewer));
487:   PetscFunctionReturn(PETSC_SUCCESS);
488: }

490: /*@
491:   PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing.

493:   Not Collective

495:   Input Parameter:
496: . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`

498:   Level: intermediate

500:   Note:
501:   See `PetscViewerHDF5PushTimestepping()` for details.

503: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
504: @*/
505: PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer)
506: {
507:   PetscFunctionBegin;
509:   PetscUseMethod(viewer, "PetscViewerHDF5PopTimestepping_C", (PetscViewer), (viewer));
510:   PetscFunctionReturn(PETSC_SUCCESS);
511: }

513: /*@
514:   PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently.

516:   Not Collective

518:   Input Parameter:
519: . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`

521:   Output Parameter:
522: . flg - is timestepping active?

524:   Level: intermediate

526:   Note:
527:   See `PetscViewerHDF5PushTimestepping()` for details.

529: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
530: @*/
531: PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg)
532: {
533:   PetscFunctionBegin;
535:   PetscUseMethod(viewer, "PetscViewerHDF5IsTimestepping_C", (PetscViewer, PetscBool *), (viewer, flg));
536:   PetscFunctionReturn(PETSC_SUCCESS);
537: }

539: /*@
540:   PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time.

542:   Not Collective

544:   Input Parameter:
545: . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`

547:   Level: intermediate

549:   Note:
550:   This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.

552: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5GetTimestep()`
553: @*/
554: PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
555: {
556:   PetscFunctionBegin;
558:   PetscUseMethod(viewer, "PetscViewerHDF5IncrementTimestep_C", (PetscViewer), (viewer));
559:   PetscFunctionReturn(PETSC_SUCCESS);
560: }

562: /*@
563:   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time.

565:   Logically Collective

567:   Input Parameters:
568: + viewer   - the `PetscViewer` of type `PETSCVIEWERHDF5`
569: - timestep - The timestep

571:   Level: intermediate

573:   Note:
574:   This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.

576: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
577: @*/
578: PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
579: {
580:   PetscFunctionBegin;
583:   PetscUseMethod(viewer, "PetscViewerHDF5SetTimestep_C", (PetscViewer, PetscInt), (viewer, timestep));
584:   PetscFunctionReturn(PETSC_SUCCESS);
585: }

587: /*@
588:   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.

590:   Not Collective

592:   Input Parameter:
593: . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`

595:   Output Parameter:
596: . timestep - The timestep

598:   Level: intermediate

600:   Note:
601:   This can be called only if the viewer is in the timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.

603: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5SetTimestep()`
604: @*/
605: PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
606: {
607:   PetscFunctionBegin;
609:   PetscAssertPointer(timestep, 2);
610:   PetscUseMethod(viewer, "PetscViewerHDF5GetTimestep_C", (PetscViewer, PetscInt *), (viewer, timestep));
611:   PetscFunctionReturn(PETSC_SUCCESS);
612: }

614: /*@C
615:   PetscViewerHDF5WriteAttribute - Write an attribute

617:   Collective

619:   Input Parameters:
620: + viewer   - The `PETSCVIEWERHDF5` viewer
621: . parent   - The parent dataset/group name
622: . name     - The attribute name
623: . datatype - The attribute type
624: - value    - The attribute value

626:   Level: advanced

628:   Note:
629:   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group.
630:   `NULL` means the current pushed group.

632: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5HasAttribute()`,
633:           `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
634: @*/
635: PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
636: {
637:   PetscFunctionBegin;
639:   if (parent) PetscAssertPointer(parent, 2);
640:   PetscAssertPointer(name, 3);
642:   PetscAssertPointer(value, 5);
643:   PetscUseMethod(viewer, "PetscViewerHDF5WriteAttribute_C", (PetscViewer, const char[], const char[], PetscDataType, const void *), (viewer, parent, name, datatype, value));
644:   PetscFunctionReturn(PETSC_SUCCESS);
645: }

647: /*@C
648:   PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given `PetscObject` by name

650:   Collective

652:   Input Parameters:
653: + viewer   - The `PETSCVIEWERHDF5` viewer
654: . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
655: . name     - The attribute name
656: . datatype - The attribute type
657: - value    - The attribute value

659:   Level: advanced

661:   Note:
662:   This fails if the path current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
663:   You might want to check first if it does using `PetscViewerHDF5HasObject()`.

665: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`,
666:           `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
667: @*/
668: PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
669: {
670:   PetscFunctionBegin;
673:   PetscAssertPointer(name, 3);
675:   PetscAssertPointer(value, 5);
676:   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
677:   PetscCall(PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value));
678:   PetscFunctionReturn(PETSC_SUCCESS);
679: }

681: /*@C
682:   PetscViewerHDF5ReadAttribute - Read an attribute

684:   Collective

686:   Input Parameters:
687: + viewer       - The `PETSCVIEWERHDF5` viewer
688: . parent       - The parent dataset/group name
689: . name         - The attribute name
690: . datatype     - The attribute type
691: - defaultValue - The pointer to the default value

693:   Output Parameter:
694: . value - The pointer to the read HDF5 attribute value

696:   Level: advanced

698:   Notes:
699:   If defaultValue is `NULL` and the attribute is not found, an error occurs.

701:   If defaultValue is not `NULL` and the attribute is not found, `defaultValue` is copied to value.

703:   The pointers `defaultValue` and `value` can be the same; for instance
704: .vb
705:   flg = PETSC_FALSE;
706:   PetscCall(`PetscViewerHDF5ReadAttribute`(viewer,name,"attr",PETSC_BOOL,&flg,&flg));
707: .ve
708:   is valid, but make sure the default value is initialized.

710:   If the datatype is `PETSC_STRING`, the output string is newly allocated so one must `PetscFree()` it when no longer needed.

712:   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group. `NULL` means the current pushed group.

714: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
715: @*/
716: PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value)
717: {
718:   PetscFunctionBegin;
720:   if (parent) PetscAssertPointer(parent, 2);
721:   PetscAssertPointer(name, 3);
723:   if (defaultValue) PetscAssertPointer(defaultValue, 5);
724:   PetscAssertPointer(value, 6);
725:   PetscUseMethod(viewer, "PetscViewerHDF5ReadAttribute_C", (PetscViewer, const char[], const char[], PetscDataType, const void *, void *), (viewer, parent, name, datatype, defaultValue, value));
726:   PetscFunctionReturn(PETSC_SUCCESS);
727: }

729: /*@C
730:   PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given `PetscObject` by name

732:   Collective

734:   Input Parameters:
735: + viewer       - The `PETSCVIEWERHDF5` viewer
736: . obj          - The object whose name is used to lookup the parent dataset, relative to the current group.
737: . name         - The attribute name
738: . datatype     - The attribute type
739: - defaultValue - The default attribute value

741:   Output Parameter:
742: . value - The attribute value

744:   Level: advanced

746:   Note:
747:   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
748:   You might want to check first if it does using `PetscViewerHDF5HasObject()`.

750: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`,
751:           `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
752: @*/
753: PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value)
754: {
755:   PetscFunctionBegin;
758:   PetscAssertPointer(name, 3);
759:   PetscAssertPointer(value, 6);
760:   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
761:   PetscCall(PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value));
762:   PetscFunctionReturn(PETSC_SUCCESS);
763: }

765: /*@C
766:   PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file

768:   Collective

770:   Input Parameters:
771: + viewer - The `PETSCVIEWERHDF5` viewer
772: - path   - (Optional) The path relative to the pushed group

774:   Output Parameter:
775: . has - Flag for group existence

777:   Level: advanced

779:   Notes:
780:   If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group.
781:   `NULL` or empty path means the current pushed group.

783:   If path exists but is not a group, `PETSC_FALSE` is returned.

785: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
786: @*/
787: PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has)
788: {
789:   PetscFunctionBegin;
791:   if (path) PetscAssertPointer(path, 2);
792:   PetscAssertPointer(has, 3);
793:   PetscUseMethod(viewer, "PetscViewerHDF5HasGroup_C", (PetscViewer, const char[], PetscBool *), (viewer, path, has));
794:   PetscFunctionReturn(PETSC_SUCCESS);
795: }

797: /*@C
798:   PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file

800:   Collective

802:   Input Parameters:
803: + viewer - The `PETSCVIEWERHDF5` viewer
804: - path   - The dataset path

806:   Output Parameter:
807: . has - Flag whether dataset exists

809:   Level: advanced

811:   Notes:
812:   If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group.

814:   If `path` is `NULL` or empty, has is set to `PETSC_FALSE`.

816:   If `path` exists but is not a dataset, has is set to `PETSC_FALSE` as well.

818: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasGroup()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
819: @*/
820: PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has)
821: {
822:   PetscFunctionBegin;
824:   if (path) PetscAssertPointer(path, 2);
825:   PetscAssertPointer(has, 3);
826:   PetscUseMethod(viewer, "PetscViewerHDF5HasDataset_C", (PetscViewer, const char[], PetscBool *), (viewer, path, has));
827:   PetscFunctionReturn(PETSC_SUCCESS);
828: }

830: /*@
831:   PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group

833:   Collective

835:   Input Parameters:
836: + viewer - The `PETSCVIEWERHDF5` viewer
837: - obj    - The named object

839:   Output Parameter:
840: . has - Flag for dataset existence

842:   Level: advanced

844:   Notes:
845:   If the object is unnamed, an error occurs.

847:   If the path current_group/object_name exists but is not a dataset, has is set to `PETSC_FALSE` as well.

849: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
850: @*/
851: PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
852: {
853:   size_t len;

855:   PetscFunctionBegin;
858:   PetscAssertPointer(has, 3);
859:   PetscCall(PetscStrlen(obj->name, &len));
860:   PetscCheck(len, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
861:   PetscCall(PetscViewerHDF5HasDataset(viewer, obj->name, has));
862:   PetscFunctionReturn(PETSC_SUCCESS);
863: }

865: /*@C
866:   PetscViewerHDF5HasAttribute - Check whether an attribute exists

868:   Collective

870:   Input Parameters:
871: + viewer - The `PETSCVIEWERHDF5` viewer
872: . parent - The parent dataset/group name
873: - name   - The attribute name

875:   Output Parameter:
876: . has - Flag for attribute existence

878:   Level: advanced

880:   Note:
881:   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group. `NULL` means the current pushed group.

883: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
884: @*/
885: PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
886: {
887:   PetscFunctionBegin;
889:   if (parent) PetscAssertPointer(parent, 2);
890:   PetscAssertPointer(name, 3);
891:   PetscAssertPointer(has, 4);
892:   PetscUseMethod(viewer, "PetscViewerHDF5HasAttribute_C", (PetscViewer, const char[], const char[], PetscBool *), (viewer, parent, name, has));
893:   PetscFunctionReturn(PETSC_SUCCESS);
894: }

896: /*@C
897:   PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given `PetscObject` by name

899:   Collective

901:   Input Parameters:
902: + viewer - The `PETSCVIEWERHDF5` viewer
903: . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
904: - name   - The attribute name

906:   Output Parameter:
907: . has - Flag for attribute existence

909:   Level: advanced

911:   Note:
912:   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
913:   You might want to check first if it does using `PetscViewerHDF5HasObject()`.

915: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
916: @*/
917: PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
918: {
919:   PetscFunctionBegin;
922:   PetscAssertPointer(name, 3);
923:   PetscAssertPointer(has, 4);
924:   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
925:   PetscCall(PetscViewerHDF5HasAttribute(viewer, obj->name, name, has));
926:   PetscFunctionReturn(PETSC_SUCCESS);
927: }

929: /*
930:   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
931:   is attached to a communicator, in this case the attribute is a PetscViewer.
932: */
933: PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;

935: /*@C
936:   PETSC_VIEWER_HDF5_ - Creates an `PETSCVIEWERHDF5` `PetscViewer` shared by all processors in a communicator.

938:   Collective

940:   Input Parameter:
941: . comm - the MPI communicator to share the `PETSCVIEWERHDF5` `PetscViewer`

943:   Options Database Key:
944: . -viewer_hdf5_filename name - name of the HDF5 file

946:   Environmental variable:
947: . `PETSC_VIEWER_HDF5_FILENAME` - name of the HDF5 file

949:   Level: intermediate

951:   Note:
952:   Unlike almost all other PETSc routines, `PETSC_VIEWER_HDF5_()` does not return
953:   an error code.  The HDF5 `PetscViewer` is usually used in the form
954: .vb
955:   XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
956: .ve

958: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerCreate()`, `PetscViewerDestroy()`
959: @*/
960: PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
961: {
962:   PetscErrorCode ierr;
963:   PetscMPIInt    mpi_ierr;
964:   PetscBool      flg;
965:   PetscMPIInt    iflg;
966:   PetscViewer    viewer;
967:   char           fname[PETSC_MAX_PATH_LEN];
968:   MPI_Comm       ncomm;

970:   PetscFunctionBegin;
971:   ierr = PetscCommDuplicate(comm, &ncomm, NULL);
972:   if (ierr) {
973:     ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
974:     PetscFunctionReturn(NULL);
975:   }
976:   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
977:     mpi_ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Viewer_HDF5_keyval, NULL);
978:     if (mpi_ierr) {
979:       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
980:       PetscFunctionReturn(NULL);
981:     }
982:   }
983:   mpi_ierr = MPI_Comm_get_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void **)&viewer, &iflg);
984:   if (mpi_ierr) {
985:     ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
986:     PetscFunctionReturn(NULL);
987:   }
988:   if (!iflg) { /* PetscViewer not yet created */
989:     ierr = PetscOptionsGetenv(ncomm, "PETSC_VIEWER_HDF5_FILENAME", fname, PETSC_MAX_PATH_LEN, &flg);
990:     if (ierr) {
991:       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
992:       PetscFunctionReturn(NULL);
993:     }
994:     if (!flg) {
995:       ierr = PetscStrncpy(fname, "output.h5", sizeof(fname));
996:       if (ierr) {
997:         ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
998:         PetscFunctionReturn(NULL);
999:       }
1000:     }
1001:     ierr = PetscViewerHDF5Open(ncomm, fname, FILE_MODE_WRITE, &viewer);
1002:     if (ierr) {
1003:       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
1004:       PetscFunctionReturn(NULL);
1005:     }
1006:     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
1007:     if (ierr) {
1008:       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
1009:       PetscFunctionReturn(NULL);
1010:     }
1011:     mpi_ierr = MPI_Comm_set_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void *)viewer);
1012:     if (mpi_ierr) {
1013:       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
1014:       PetscFunctionReturn(NULL);
1015:     }
1016:   }
1017:   ierr = PetscCommDestroy(&ncomm);
1018:   if (ierr) {
1019:     ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
1020:     PetscFunctionReturn(NULL);
1021:   }
1022:   PetscFunctionReturn(viewer);
1023: }