Actual source code: hdf5v.c

  1: #include <petsc/private/viewerhdf5impl.h>

  3: static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool *, H5O_type_t *);
  4: static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool *);
  5: static PetscErrorCode PetscViewerHDF5GetGroup_Internal(PetscViewer, const char *[]);

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

 10:   Not Collective

 12:   Input Parameters:
 13: + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
 14: - path   - (Optional) The path relative to the pushed group

 16:   Output Parameter:
 17: . abspath - The absolute HDF5 path (group)

 19:   Level: intermediate

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

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

 27: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()`
 28: @*/
 29: PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char path[], const char *abspath[])
 30: {
 31:   size_t      len;
 32:   PetscBool   relative = PETSC_FALSE;
 33:   const char *group;
 34:   char        buf[PETSC_MAX_PATH_LEN] = "";

 36:   PetscFunctionBegin;
 38:   if (path) PetscAssertPointer(path, 2);
 39:   PetscAssertPointer(abspath, 3);
 40:   PetscCall(PetscViewerHDF5GetGroup_Internal(viewer, &group));
 41:   PetscCall(PetscStrlen(path, &len));
 42:   relative = (PetscBool)(!len || path[0] != '/');
 43:   if (relative) {
 44:     PetscCall(PetscStrncpy(buf, group, sizeof(buf)));
 45:     if (!group || len) PetscCall(PetscStrlcat(buf, "/", sizeof(buf)));
 46:     PetscCall(PetscStrlcat(buf, path, sizeof(buf)));
 47:     PetscCall(PetscStrallocpy(buf, (char **)abspath));
 48:   } else {
 49:     PetscCall(PetscStrallocpy(path, (char **)abspath));
 50:   }
 51:   PetscFunctionReturn(PETSC_SUCCESS);
 52: }

 54: static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj)
 55: {
 56:   PetscBool has;

 58:   PetscFunctionBegin;
 59:   PetscCall(PetscViewerHDF5HasObject(viewer, obj, &has));
 60:   if (!has) {
 61:     const char *group;
 62:     PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group));
 63:     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", obj->name, group);
 64:   }
 65:   PetscFunctionReturn(PETSC_SUCCESS);
 66: }

 68: static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscViewer v, PetscOptionItems PetscOptionsObject)
 69: {
 70:   PetscBool         flg  = PETSC_FALSE, set;
 71:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data;

 73:   PetscFunctionBegin;
 74:   PetscOptionsHeadBegin(PetscOptionsObject, "HDF5 PetscViewer Options");
 75:   PetscCall(PetscOptionsBool("-viewer_hdf5_base_dimension2", "1d Vectors get 2 dimensions in HDF5", "PetscViewerHDF5SetBaseDimension2", hdf5->basedimension2, &hdf5->basedimension2, NULL));
 76:   PetscCall(PetscOptionsBool("-viewer_hdf5_sp_output", "Force data to be written in single precision", "PetscViewerHDF5SetSPOutput", hdf5->spoutput, &hdf5->spoutput, NULL));
 77:   PetscCall(PetscOptionsBool("-viewer_hdf5_collective", "Enable collective transfer mode", "PetscViewerHDF5SetCollective", flg, &flg, &set));
 78:   if (set) PetscCall(PetscViewerHDF5SetCollective(v, flg));
 79:   flg = PETSC_FALSE;
 80:   PetscCall(PetscOptionsBool("-viewer_hdf5_default_timestepping", "Set default timestepping state", "PetscViewerHDF5SetDefaultTimestepping", flg, &flg, &set));
 81:   if (set) PetscCall(PetscViewerHDF5SetDefaultTimestepping(v, flg));
 82:   PetscCall(PetscOptionsBool("-viewer_hdf5_compress", "Enable compression", "PetscViewerHDF5SetCompress", hdf5->compress, &hdf5->compress, NULL));
 83:   PetscOptionsHeadEnd();
 84:   PetscFunctionReturn(PETSC_SUCCESS);
 85: }

 87: static PetscErrorCode PetscViewerView_HDF5(PetscViewer v, PetscViewer viewer)
 88: {
 89:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data;
 90:   PetscBool         flg;

 92:   PetscFunctionBegin;
 93:   if (hdf5->filename) PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", hdf5->filename));
 94:   PetscCall(PetscViewerASCIIPrintf(viewer, "Vectors with blocksize 1 saved as 2D datasets: %s\n", PetscBools[hdf5->basedimension2]));
 95:   PetscCall(PetscViewerASCIIPrintf(viewer, "Enforce single precision storage: %s\n", PetscBools[hdf5->spoutput]));
 96:   PetscCall(PetscViewerHDF5GetCollective(v, &flg));
 97:   PetscCall(PetscViewerASCIIPrintf(viewer, "MPI-IO transfer mode: %s\n", flg ? "collective" : "independent"));
 98:   PetscCall(PetscViewerASCIIPrintf(viewer, "Default timestepping: %s\n", PetscBools[hdf5->defTimestepping]));
 99:   PetscCall(PetscViewerASCIIPrintf(viewer, "Compression: %s\n", PetscBools[hdf5->compress]));
100:   PetscFunctionReturn(PETSC_SUCCESS);
101: }

103: static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
104: {
105:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

107:   PetscFunctionBegin;
108:   PetscCall(PetscFree(hdf5->filename));
109:   if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id));
110:   PetscFunctionReturn(PETSC_SUCCESS);
111: }

113: static PetscErrorCode PetscViewerFlush_HDF5(PetscViewer viewer)
114: {
115:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

117:   PetscFunctionBegin;
118:   if (hdf5->file_id) PetscCallHDF5(H5Fflush, (hdf5->file_id, H5F_SCOPE_LOCAL));
119:   PetscFunctionReturn(PETSC_SUCCESS);
120: }

122: static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
123: {
124:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

126:   PetscFunctionBegin;
127:   PetscCallHDF5(H5Pclose, (hdf5->dxpl_id));
128:   PetscCall(PetscViewerFileClose_HDF5(viewer));
129:   while (hdf5->groups) {
130:     PetscViewerHDF5GroupList *tmp = hdf5->groups->next;

132:     PetscCall(PetscFree(hdf5->groups->name));
133:     PetscCall(PetscFree(hdf5->groups));
134:     hdf5->groups = tmp;
135:   }
136:   PetscCall(PetscFree(hdf5));
137:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetName_C", NULL));
138:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetName_C", NULL));
139:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetMode_C", NULL));
140:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetMode_C", NULL));
141:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetBaseDimension2_C", NULL));
142:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetSPOutput_C", NULL));
143:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetCollective_C", NULL));
144:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetCollective_C", NULL));
145:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetDefaultTimestepping_C", NULL));
146:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetDefaultTimestepping_C", NULL));
147:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetCompress_C", NULL));
148:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetCompress_C", NULL));
149:   PetscFunctionReturn(PETSC_SUCCESS);
150: }

152: static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
153: {
154:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

156:   PetscFunctionBegin;
157:   hdf5->btype = type;
158:   PetscFunctionReturn(PETSC_SUCCESS);
159: }

161: static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
162: {
163:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

165:   PetscFunctionBegin;
166:   *type = hdf5->btype;
167:   PetscFunctionReturn(PETSC_SUCCESS);
168: }

170: static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
171: {
172:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

174:   PetscFunctionBegin;
175:   hdf5->basedimension2 = flg;
176:   PetscFunctionReturn(PETSC_SUCCESS);
177: }

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

183:   Logically Collective

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

189:   Options Database Key:
190: . -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

192:   Level: intermediate

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

198: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
199: @*/
200: PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer, PetscBool flg)
201: {
202:   PetscFunctionBegin;
204:   PetscTryMethod(viewer, "PetscViewerHDF5SetBaseDimension2_C", (PetscViewer, PetscBool), (viewer, flg));
205:   PetscFunctionReturn(PETSC_SUCCESS);
206: }

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

212:   Logically Collective

214:   Input Parameter:
215: . viewer - the `PetscViewer`, must be `PETSCVIEWERHDF5`

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

220:   Level: intermediate

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

226: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
227: @*/
228: PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer, PetscBool *flg)
229: {
230:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

232:   PetscFunctionBegin;
234:   *flg = hdf5->basedimension2;
235:   PetscFunctionReturn(PETSC_SUCCESS);
236: }

238: static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
239: {
240:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

242:   PetscFunctionBegin;
243:   hdf5->spoutput = flg;
244:   PetscFunctionReturn(PETSC_SUCCESS);
245: }

247: /*@
248:   PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
249:   compiled with double precision `PetscReal`.

251:   Logically Collective

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

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

260:   Level: intermediate

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

266: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
267:           `PetscReal`, `PetscViewerHDF5GetSPOutput()`
268: @*/
269: PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer, PetscBool flg)
270: {
271:   PetscFunctionBegin;
273:   PetscTryMethod(viewer, "PetscViewerHDF5SetSPOutput_C", (PetscViewer, PetscBool), (viewer, flg));
274:   PetscFunctionReturn(PETSC_SUCCESS);
275: }

277: /*@
278:   PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
279:   compiled with double precision `PetscReal`.

281:   Logically Collective

283:   Input Parameter:
284: . viewer - the PetscViewer, must be of type `PETSCVIEWERHDF5`

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

289:   Level: intermediate

291: .seealso: [](sec_viewers), `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
292:           `PetscReal`, `PetscViewerHDF5SetSPOutput()`
293: @*/
294: PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer, PetscBool *flg)
295: {
296:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

298:   PetscFunctionBegin;
300:   *flg = hdf5->spoutput;
301:   PetscFunctionReturn(PETSC_SUCCESS);
302: }

304: static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg)
305: {
306:   PetscFunctionBegin;
307:   /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions
308:      - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */
309: #if H5_VERSION_GE(1, 10, 3) && defined(H5_HAVE_PARALLEL)
310:   {
311:     PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
312:     PetscCallHDF5(H5Pset_dxpl_mpio, (hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT));
313:   }
314: #else
315:   if (flg) PetscCall(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"));
316: #endif
317:   PetscFunctionReturn(PETSC_SUCCESS);
318: }

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

323:   Logically Collective; flg must contain common value

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

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

332:   Level: intermediate

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

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

343: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5GetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()`
344: @*/
345: PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer, PetscBool flg)
346: {
347:   PetscFunctionBegin;
350:   PetscTryMethod(viewer, "PetscViewerHDF5SetCollective_C", (PetscViewer, PetscBool), (viewer, flg));
351:   PetscFunctionReturn(PETSC_SUCCESS);
352: }

354: static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg)
355: {
356: #if defined(H5_HAVE_PARALLEL)
357:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
358:   H5FD_mpio_xfer_t  mode;
359: #endif

361:   PetscFunctionBegin;
362: #if !defined(H5_HAVE_PARALLEL)
363:   *flg = PETSC_FALSE;
364: #else
365:   PetscCallHDF5(H5Pget_dxpl_mpio, (hdf5->dxpl_id, &mode));
366:   *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE;
367: #endif
368:   PetscFunctionReturn(PETSC_SUCCESS);
369: }

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

374:   Not Collective

376:   Input Parameter:
377: . viewer - the `PETSCVIEWERHDF5` `PetscViewer`

379:   Output Parameter:
380: . flg - the flag

382:   Level: intermediate

384:   Note:
385:   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.
386:   For more details, see `PetscViewerHDF5SetCollective()`.

388: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5SetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()`
389: @*/
390: PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer, PetscBool *flg)
391: {
392:   PetscFunctionBegin;
394:   PetscAssertPointer(flg, 2);

396:   PetscUseMethod(viewer, "PetscViewerHDF5GetCollective_C", (PetscViewer, PetscBool *), (viewer, flg));
397:   PetscFunctionReturn(PETSC_SUCCESS);
398: }

400: static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
401: {
402:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
403:   hid_t             plist_id;

405:   PetscFunctionBegin;
406:   if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id));
407:   if (hdf5->filename) PetscCall(PetscFree(hdf5->filename));
408:   PetscCall(PetscStrallocpy(name, &hdf5->filename));
409:   /* Set up file access property list with parallel I/O access */
410:   PetscCallHDF5Return(plist_id, H5Pcreate, (H5P_FILE_ACCESS));
411: #if defined(H5_HAVE_PARALLEL)
412:   PetscCallHDF5(H5Pset_fapl_mpio, (plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL));
413: #endif
414:   /* Create or open the file collectively */
415:   switch (hdf5->btype) {
416:   case FILE_MODE_READ:
417:     if (PetscDefined(USE_DEBUG)) {
418:       PetscMPIInt rank;
419:       PetscBool   flg;

421:       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
422:       if (rank == 0) {
423:         PetscCall(PetscTestFile(hdf5->filename, 'r', &flg));
424:         PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "File %s requested for reading does not exist", hdf5->filename);
425:       }
426:       PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer)));
427:     }
428:     PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDONLY, plist_id));
429:     break;
430:   case FILE_MODE_APPEND:
431:   case FILE_MODE_UPDATE: {
432:     PetscBool flg;
433:     PetscCall(PetscTestFile(hdf5->filename, 'r', &flg));
434:     if (flg) PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDWR, plist_id));
435:     else PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id));
436:     break;
437:   }
438:   case FILE_MODE_WRITE:
439:     PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
440:     break;
441:   case FILE_MODE_UNDEFINED:
442:     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
443:   default:
444:     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[hdf5->btype]);
445:   }
446:   PetscCheck(hdf5->file_id >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
447:   PetscCallHDF5(H5Pclose, (plist_id));
448:   PetscCall(PetscViewerHDF5ResetAttachedDMPlexStorageVersion(viewer));
449:   PetscFunctionReturn(PETSC_SUCCESS);
450: }

452: static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer, const char **name)
453: {
454:   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5 *)viewer->data;

456:   PetscFunctionBegin;
457:   *name = vhdf5->filename;
458:   PetscFunctionReturn(PETSC_SUCCESS);
459: }

461: static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
462: {
463:   /*
464:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
465:   */

467:   PetscFunctionBegin;
468:   PetscFunctionReturn(PETSC_SUCCESS);
469: }

471: static PetscErrorCode PetscViewerHDF5SetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool flg)
472: {
473:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

475:   PetscFunctionBegin;
476:   hdf5->defTimestepping = flg;
477:   PetscFunctionReturn(PETSC_SUCCESS);
478: }

480: static PetscErrorCode PetscViewerHDF5GetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool *flg)
481: {
482:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

484:   PetscFunctionBegin;
485:   *flg = hdf5->defTimestepping;
486:   PetscFunctionReturn(PETSC_SUCCESS);
487: }

489: /*@
490:   PetscViewerHDF5SetDefaultTimestepping - Set the flag for default timestepping

492:   Logically Collective

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

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

501:   Level: intermediate

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

506: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5GetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()`
507: @*/
508: PetscErrorCode PetscViewerHDF5SetDefaultTimestepping(PetscViewer viewer, PetscBool flg)
509: {
510:   PetscFunctionBegin;
512:   PetscTryMethod(viewer, "PetscViewerHDF5SetDefaultTimestepping_C", (PetscViewer, PetscBool), (viewer, flg));
513:   PetscFunctionReturn(PETSC_SUCCESS);
514: }

516: /*@
517:   PetscViewerHDF5GetDefaultTimestepping - Get the flag for default timestepping

519:   Not Collective

521:   Input Parameter:
522: . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`

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

527:   Level: intermediate

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

539: static PetscErrorCode PetscViewerHDF5SetCompress_HDF5(PetscViewer viewer, PetscBool flg)
540: {
541:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

543:   PetscFunctionBegin;
544:   hdf5->compress = flg;
545:   PetscFunctionReturn(PETSC_SUCCESS);
546: }

548: static PetscErrorCode PetscViewerHDF5GetCompress_HDF5(PetscViewer viewer, PetscBool *flg)
549: {
550:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

552:   PetscFunctionBegin;
553:   *flg = hdf5->compress;
554:   PetscFunctionReturn(PETSC_SUCCESS);
555: }

557: /*@
558:   PetscViewerHDF5SetCompress - Set the flag for compression

560:   Logically Collective

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

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

569:   Level: intermediate

571: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5GetCompress()`
572: @*/
573: PetscErrorCode PetscViewerHDF5SetCompress(PetscViewer viewer, PetscBool flg)
574: {
575:   PetscFunctionBegin;
577:   PetscTryMethod(viewer, "PetscViewerHDF5SetCompress_C", (PetscViewer, PetscBool), (viewer, flg));
578:   PetscFunctionReturn(PETSC_SUCCESS);
579: }

581: /*@
582:   PetscViewerHDF5GetCompress - Get the flag for compression

584:   Not Collective

586:   Input Parameter:
587: . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`

589:   Output Parameter:
590: . flg - if `PETSC_TRUE` we will turn on compression

592:   Level: intermediate

594: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5SetCompress()`
595: @*/
596: PetscErrorCode PetscViewerHDF5GetCompress(PetscViewer viewer, PetscBool *flg)
597: {
598:   PetscFunctionBegin;
600:   PetscUseMethod(viewer, "PetscViewerHDF5GetCompress_C", (PetscViewer, PetscBool *), (viewer, flg));
601:   PetscFunctionReturn(PETSC_SUCCESS);
602: }

604: /*MC
605:    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file

607:   Level: beginner

609: .seealso: [](sec_viewers), `PetscViewerHDF5Open()`, `PetscViewerStringSPrintf()`, `PetscViewerSocketOpen()`, `PetscViewerDrawOpen()`, `PETSCVIEWERSOCKET`,
610:           `PetscViewerCreate()`, `PetscViewerASCIIOpen()`, `PetscViewerBinaryOpen()`, `PETSCVIEWERBINARY`, `PETSCVIEWERDRAW`, `PETSCVIEWERSTRING`,
611:           `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`,
612:           `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`
613: M*/

615: PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
616: {
617:   PetscViewer_HDF5 *hdf5;

619:   PetscFunctionBegin;
620: #if !defined(H5_HAVE_PARALLEL)
621:   {
622:     PetscMPIInt size;
623:     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)v), &size));
624:     PetscCheck(size <= 1, PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "Cannot use parallel HDF5 viewer since the given HDF5 does not support parallel I/O (H5_HAVE_PARALLEL is unset)");
625:   }
626: #endif

628:   PetscCall(PetscNew(&hdf5));

630:   v->data                = (void *)hdf5;
631:   v->ops->destroy        = PetscViewerDestroy_HDF5;
632:   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
633:   v->ops->setup          = PetscViewerSetUp_HDF5;
634:   v->ops->view           = PetscViewerView_HDF5;
635:   v->ops->flush          = PetscViewerFlush_HDF5;
636:   hdf5->btype            = FILE_MODE_UNDEFINED;
637:   hdf5->filename         = NULL;
638:   hdf5->timestep         = -1;
639:   hdf5->groups           = NULL;

641:   PetscCallHDF5Return(hdf5->dxpl_id, H5Pcreate, (H5P_DATASET_XFER));

643:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_HDF5));
644:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_HDF5));
645:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_HDF5));
646:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_HDF5));
647:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetBaseDimension2_C", PetscViewerHDF5SetBaseDimension2_HDF5));
648:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetSPOutput_C", PetscViewerHDF5SetSPOutput_HDF5));
649:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetCollective_C", PetscViewerHDF5SetCollective_HDF5));
650:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetCollective_C", PetscViewerHDF5GetCollective_HDF5));
651:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetDefaultTimestepping_C", PetscViewerHDF5GetDefaultTimestepping_HDF5));
652:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetDefaultTimestepping_C", PetscViewerHDF5SetDefaultTimestepping_HDF5));
653:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetCompress_C", PetscViewerHDF5GetCompress_HDF5));
654:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetCompress_C", PetscViewerHDF5SetCompress_HDF5));
655:   PetscFunctionReturn(PETSC_SUCCESS);
656: }

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

661:   Collective

663:   Input Parameters:
664: + comm - MPI communicator
665: . name - name of file
666: - type - type of file

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

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

675:   Level: beginner

677:   Notes:
678:   Reading is always available, regardless of the mode. Available modes are
679: .vb
680:   FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY]
681:   FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC]
682:   FILE_MODE_APPEND - if file exists, keep existing contents [H5Fopen() with H5F_ACC_RDWR], else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_EXCL]
683:   FILE_MODE_UPDATE - same as FILE_MODE_APPEND
684: .ve

686:   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.

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

690: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, `PetscViewerHDF5SetBaseDimension2()`,
691:           `PetscViewerHDF5SetSPOutput()`, `PetscViewerHDF5GetBaseDimension2()`, `VecView()`, `MatView()`, `VecLoad()`,
692:           `MatLoad()`, `PetscFileMode`, `PetscViewer`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()`
693: @*/
694: PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
695: {
696:   PetscFunctionBegin;
697:   PetscCall(PetscViewerCreate(comm, hdf5v));
698:   PetscCall(PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5));
699:   PetscCall(PetscViewerFileSetMode(*hdf5v, type));
700:   PetscCall(PetscViewerFileSetName(*hdf5v, name));
701:   PetscCall(PetscViewerSetFromOptions(*hdf5v));
702:   PetscFunctionReturn(PETSC_SUCCESS);
703: }

705: /*@C
706:   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls

708:   Not Collective

710:   Input Parameter:
711: . viewer - the `PetscViewer`

713:   Output Parameter:
714: . file_id - The file id

716:   Level: intermediate

718: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`
719: @*/
720: PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
721: {
722:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

724:   PetscFunctionBegin;
726:   if (file_id) *file_id = hdf5->file_id;
727:   PetscFunctionReturn(PETSC_SUCCESS);
728: }

730: /*@C
731:   PetscViewerHDF5PushGroup - Set the current HDF5 group for output

733:   Not Collective

735:   Input Parameters:
736: + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
737: - name   - The group name

739:   Level: intermediate

741:   Notes:
742:   This is designed to mnemonically resemble the Unix cd command.
743: .vb
744:   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.
745:   `NULL`, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/".
746:   "." means the current group is pushed again.
747: .ve

749:   Example:
750:   Suppose the current group is "/a".
751: .vb
752:   If name is `NULL`, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/".
753:   If name is ".", then the new group will be "/a".
754:   If name is "b", then the new group will be "/a/b".
755:   If name is "/b", then the new group will be "/b".
756: .ve

758:   Developer Notes:
759:   The root group "/" is internally stored as `NULL`.

761: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()`
762: @*/
763: PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[])
764: {
765:   PetscViewer_HDF5         *hdf5 = (PetscViewer_HDF5 *)viewer->data;
766:   PetscViewerHDF5GroupList *groupNode;
767:   size_t                    i, len;
768:   char                      buf[PETSC_MAX_PATH_LEN];
769:   const char               *gname;

771:   PetscFunctionBegin;
773:   if (name) PetscAssertPointer(name, 2);
774:   PetscCall(PetscStrlen(name, &len));
775:   gname = NULL;
776:   if (len) {
777:     if (len == 1 && name[0] == '.') {
778:       /* use current name */
779:       gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL;
780:     } else if (name[0] == '/') {
781:       /* absolute */
782:       for (i = 1; i < len; i++) {
783:         if (name[i] != '/') {
784:           gname = name;
785:           break;
786:         }
787:       }
788:     } else {
789:       /* relative */
790:       const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : "";
791:       PetscCall(PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name));
792:       gname = buf;
793:     }
794:   }
795:   PetscCall(PetscNew(&groupNode));
796:   PetscCall(PetscStrallocpy(gname, (char **)&groupNode->name));
797:   groupNode->next = hdf5->groups;
798:   hdf5->groups    = groupNode;
799:   PetscFunctionReturn(PETSC_SUCCESS);
800: }

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

805:   Not Collective

807:   Input Parameter:
808: . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`

810:   Level: intermediate

812: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()`
813: @*/
814: PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer)
815: {
816:   PetscViewer_HDF5         *hdf5 = (PetscViewer_HDF5 *)viewer->data;
817:   PetscViewerHDF5GroupList *groupNode;

819:   PetscFunctionBegin;
821:   PetscCheck(hdf5->groups, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
822:   groupNode    = hdf5->groups;
823:   hdf5->groups = hdf5->groups->next;
824:   PetscCall(PetscFree(groupNode->name));
825:   PetscCall(PetscFree(groupNode));
826:   PetscFunctionReturn(PETSC_SUCCESS);
827: }

829: static PetscErrorCode PetscViewerHDF5GetGroup_Internal(PetscViewer viewer, const char *name[])
830: {
831:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

833:   PetscFunctionBegin;
835:   PetscAssertPointer(name, 2);
836:   if (hdf5->groups) *name = hdf5->groups->name;
837:   else *name = NULL;
838:   PetscFunctionReturn(PETSC_SUCCESS);
839: }

841: /*@C
842:   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by `PetscViewerHDF5GetGroup()`,
843:   and return this group's ID and file ID.
844:   If `PetscViewerHDF5GetGroup()` yields NULL, then group ID is file ID.

846:   Not Collective

848:   Input Parameters:
849: + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
850: - path   - (Optional) The path relative to the pushed group

852:   Output Parameters:
853: + fileId  - The HDF5 file ID
854: - groupId - The HDF5 group ID

856:   Level: intermediate

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

862:   If the viewer is writable, the group is created if it doesn't exist yet.

864: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5WriteGroup()`
865: @*/
866: PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, const char path[], hid_t *fileId, hid_t *groupId)
867: {
868:   hid_t       file_id;
869:   H5O_type_t  type;
870:   const char *fileName  = NULL;
871:   const char *groupName = NULL;
872:   PetscBool   writable, has;

874:   PetscFunctionBegin;
875:   PetscCall(PetscViewerWritable(viewer, &writable));
876:   PetscCall(PetscViewerHDF5GetFileId(viewer, &file_id));
877:   PetscCall(PetscViewerFileGetName(viewer, &fileName));
878:   PetscCall(PetscViewerHDF5GetGroup(viewer, path, &groupName));
879:   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type));
880:   if (!has) {
881:     PetscCheck(writable, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Group %s does not exist and file %s is not open for writing", groupName, fileName);
882:     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName);
883:   }
884:   PetscCheck(type == H5O_TYPE_GROUP, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s in file %s resolves to something which is not a group", groupName, fileName);
885:   PetscCallHDF5Return(*groupId, H5Gopen2, (file_id, groupName, H5P_DEFAULT));
886:   PetscCall(PetscFree(groupName));
887:   *fileId = file_id;
888:   PetscFunctionReturn(PETSC_SUCCESS);
889: }

891: /*@C
892:   PetscViewerHDF5WriteGroup - Ensure the HDF5 group exists in the HDF5 file

894:   Not Collective

896:   Input Parameters:
897: + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
898: - path   - (Optional) The path relative to the pushed group

900:   Level: intermediate

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

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

908: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
909: @*/
910: PetscErrorCode PetscViewerHDF5WriteGroup(PetscViewer viewer, const char path[])
911: {
912:   hid_t fileId, groupId;

914:   PetscFunctionBegin;
915:   PetscCall(PetscViewerHDF5OpenGroup(viewer, path, &fileId, &groupId)); // make sure group is actually created
916:   PetscCallHDF5(H5Gclose, (groupId));
917:   PetscFunctionReturn(PETSC_SUCCESS);
918: }

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

923:   Not Collective

925:   Input Parameter:
926: . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`

928:   Level: intermediate

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

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

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

943: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
944: @*/
945: PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer)
946: {
947:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

949:   PetscFunctionBegin;
951:   PetscCheck(!hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed");
952:   hdf5->timestepping = PETSC_TRUE;
953:   if (hdf5->timestep < 0) hdf5->timestep = 0;
954:   PetscFunctionReturn(PETSC_SUCCESS);
955: }

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

960:   Not Collective

962:   Input Parameter:
963: . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`

965:   Level: intermediate

967:   Note:
968:   See `PetscViewerHDF5PushTimestepping()` for details.

970: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
971: @*/
972: PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer)
973: {
974:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

976:   PetscFunctionBegin;
978:   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
979:   hdf5->timestepping = PETSC_FALSE;
980:   PetscFunctionReturn(PETSC_SUCCESS);
981: }

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

986:   Not Collective

988:   Input Parameter:
989: . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`

991:   Output Parameter:
992: . flg - is timestepping active?

994:   Level: intermediate

996:   Note:
997:   See `PetscViewerHDF5PushTimestepping()` for details.

999: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
1000: @*/
1001: PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg)
1002: {
1003:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

1005:   PetscFunctionBegin;
1007:   *flg = hdf5->timestepping;
1008:   PetscFunctionReturn(PETSC_SUCCESS);
1009: }

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

1014:   Not Collective

1016:   Input Parameter:
1017: . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`

1019:   Level: intermediate

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

1024: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5GetTimestep()`
1025: @*/
1026: PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
1027: {
1028:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

1030:   PetscFunctionBegin;
1032:   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
1033:   ++hdf5->timestep;
1034:   PetscFunctionReturn(PETSC_SUCCESS);
1035: }

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

1040:   Logically Collective

1042:   Input Parameters:
1043: + viewer   - the `PetscViewer` of type `PETSCVIEWERHDF5`
1044: - timestep - The timestep

1046:   Level: intermediate

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

1051: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
1052: @*/
1053: PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
1054: {
1055:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

1057:   PetscFunctionBegin;
1060:   PetscCheck(timestep >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %" PetscInt_FMT " is negative", timestep);
1061:   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
1062:   hdf5->timestep = timestep;
1063:   PetscFunctionReturn(PETSC_SUCCESS);
1064: }

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

1069:   Not Collective

1071:   Input Parameter:
1072: . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`

1074:   Output Parameter:
1075: . timestep - The timestep

1077:   Level: intermediate

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

1082: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5SetTimestep()`
1083: @*/
1084: PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
1085: {
1086:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

1088:   PetscFunctionBegin;
1090:   PetscAssertPointer(timestep, 2);
1091:   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
1092:   *timestep = hdf5->timestep;
1093:   PetscFunctionReturn(PETSC_SUCCESS);
1094: }

1096: /*@C
1097:   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.

1099:   Not Collective

1101:   Input Parameter:
1102: . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`)

1104:   Output Parameter:
1105: . htype - the HDF5  datatype

1107:   Level: advanced

1109: .seealso: [](sec_viewers), `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()`
1110: @*/
1111: PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
1112: {
1113:   PetscFunctionBegin;
1114:   if (ptype == PETSC_INT)
1115: #if defined(PETSC_USE_64BIT_INDICES)
1116:     *htype = H5T_NATIVE_LLONG;
1117: #else
1118:     *htype = H5T_NATIVE_INT;
1119: #endif
1120:   else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE;
1121:   else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG;
1122:   else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT;
1123:   else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT;
1124:   else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT;
1125:   else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT;
1126:   else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR;
1127:   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
1128:   else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1);
1129:   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
1130:   PetscFunctionReturn(PETSC_SUCCESS);
1131: }

1133: /*@C
1134:   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name

1136:   Not Collective

1138:   Input Parameter:
1139: . htype - the HDF5 datatype (for example `H5T_NATIVE_DOUBLE`, ...)

1141:   Output Parameter:
1142: . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`)

1144:   Level: advanced

1146: .seealso: [](sec_viewers), `PetscDataType`
1147: @*/
1148: PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
1149: {
1150:   PetscFunctionBegin;
1151: #if defined(PETSC_USE_64BIT_INDICES)
1152:   if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG;
1153:   else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT;
1154: #else
1155:   if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT;
1156: #endif
1157:   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
1158:   else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG;
1159:   else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT;
1160:   else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT;
1161:   else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR;
1162:   else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR;
1163:   else if (htype == H5T_C_S1) *ptype = PETSC_STRING;
1164:   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
1165:   PetscFunctionReturn(PETSC_SUCCESS);
1166: }

1168: /*@C
1169:   PetscViewerHDF5WriteAttribute - Write an attribute

1171:   Collective

1173:   Input Parameters:
1174: + viewer   - The `PETSCVIEWERHDF5` viewer
1175: . parent   - The parent dataset/group name
1176: . name     - The attribute name
1177: . datatype - The attribute type
1178: - value    - The attribute value

1180:   Level: advanced

1182:   Note:
1183:   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.

1185: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5HasAttribute()`,
1186:           `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1187: @*/
1188: PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
1189: {
1190:   const char *parentAbsPath;
1191:   hid_t       h5, dataspace, obj, attribute, dtype;
1192:   PetscBool   has;

1194:   PetscFunctionBegin;
1196:   if (parent) PetscAssertPointer(parent, 2);
1197:   PetscAssertPointer(name, 3);
1199:   PetscAssertPointer(value, 5);
1200:   PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
1201:   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL));
1202:   PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
1203:   PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
1204:   if (datatype == PETSC_STRING) {
1205:     size_t len;
1206:     PetscCall(PetscStrlen((const char *)value, &len));
1207:     PetscCallHDF5(H5Tset_size, (dtype, len + 1));
1208:   }
1209:   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1210:   PetscCallHDF5Return(dataspace, H5Screate, (H5S_SCALAR));
1211:   PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT));
1212:   if (has) {
1213:     PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name));
1214:   } else {
1215:     PetscCallHDF5Return(attribute, H5Acreate2, (obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
1216:   }
1217:   PetscCallHDF5(H5Awrite, (attribute, dtype, value));
1218:   if (datatype == PETSC_STRING) PetscCallHDF5(H5Tclose, (dtype));
1219:   PetscCallHDF5(H5Aclose, (attribute));
1220:   PetscCallHDF5(H5Oclose, (obj));
1221:   PetscCallHDF5(H5Sclose, (dataspace));
1222:   PetscCall(PetscFree(parentAbsPath));
1223:   PetscFunctionReturn(PETSC_SUCCESS);
1224: }

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

1229:   Collective

1231:   Input Parameters:
1232: + viewer   - The `PETSCVIEWERHDF5` viewer
1233: . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1234: . name     - The attribute name
1235: . datatype - The attribute type
1236: - value    - The attribute value

1238:   Level: advanced

1240:   Note:
1241:   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).
1242:   You might want to check first if it does using `PetscViewerHDF5HasObject()`.

1244: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`,
1245:           `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1246: @*/
1247: PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
1248: {
1249:   PetscFunctionBegin;
1252:   PetscAssertPointer(name, 3);
1253:   PetscAssertPointer(value, 5);
1254:   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
1255:   PetscCall(PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value));
1256:   PetscFunctionReturn(PETSC_SUCCESS);
1257: }

1259: /*@C
1260:   PetscViewerHDF5ReadAttribute - Read an attribute

1262:   Collective

1264:   Input Parameters:
1265: + viewer       - The `PETSCVIEWERHDF5` viewer
1266: . parent       - The parent dataset/group name
1267: . name         - The attribute name
1268: . datatype     - The attribute type
1269: - defaultValue - The pointer to the default value

1271:   Output Parameter:
1272: . value - The pointer to the read HDF5 attribute value

1274:   Level: advanced

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

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

1281:   The pointers `defaultValue` and value can be the same; for instance
1282: .vb
1283:   flg = PETSC_FALSE;
1284:   PetscCall(`PetscViewerHDF5ReadAttribute`(viewer,name,"attr",PETSC_BOOL,&flg,&flg));
1285: .ve
1286:   is valid, but make sure the default value is initialized.

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

1290:   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.

1292: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1293: @*/
1294: PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value)
1295: {
1296:   const char *parentAbsPath;
1297:   hid_t       h5, obj, attribute, dtype;
1298:   PetscBool   has;

1300:   PetscFunctionBegin;
1302:   if (parent) PetscAssertPointer(parent, 2);
1303:   PetscAssertPointer(name, 3);
1304:   if (defaultValue) PetscAssertPointer(defaultValue, 5);
1305:   PetscAssertPointer(value, 6);
1306:   PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
1307:   PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
1308:   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL));
1309:   if (has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
1310:   if (!has) {
1311:     if (defaultValue) {
1312:       if (defaultValue != value) {
1313:         if (datatype == PETSC_STRING) {
1314:           PetscCall(PetscStrallocpy(*(char **)defaultValue, (char **)value));
1315:         } else {
1316:           size_t len;
1317:           PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (dtype));
1318:           PetscCall(PetscMemcpy(value, defaultValue, len));
1319:         }
1320:       }
1321:       PetscCall(PetscFree(parentAbsPath));
1322:       PetscFunctionReturn(PETSC_SUCCESS);
1323:     } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name);
1324:   }
1325:   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1326:   PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT));
1327:   PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name));
1328:   if (datatype == PETSC_STRING) {
1329:     size_t len;
1330:     hid_t  atype;
1331:     PetscCallHDF5Return(atype, H5Aget_type, (attribute));
1332:     PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (atype));
1333:     PetscCall(PetscMalloc((len + 1) * sizeof(char), value));
1334:     PetscCallHDF5(H5Tset_size, (dtype, len + 1));
1335:     PetscCallHDF5(H5Aread, (attribute, dtype, *(char **)value));
1336:   } else {
1337:     PetscCallHDF5(H5Aread, (attribute, dtype, value));
1338:   }
1339:   PetscCallHDF5(H5Aclose, (attribute));
1340:   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
1341:   PetscCallHDF5(H5Oclose, (obj));
1342:   PetscCall(PetscFree(parentAbsPath));
1343:   PetscFunctionReturn(PETSC_SUCCESS);
1344: }

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

1349:   Collective

1351:   Input Parameters:
1352: + viewer       - The `PETSCVIEWERHDF5` viewer
1353: . obj          - The object whose name is used to lookup the parent dataset, relative to the current group.
1354: . name         - The attribute name
1355: . datatype     - The attribute type
1356: - defaultValue - The default attribute value

1358:   Output Parameter:
1359: . value - The attribute value

1361:   Level: advanced

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

1367: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadAttribute()` `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1368: @*/
1369: PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value)
1370: {
1371:   PetscFunctionBegin;
1374:   PetscAssertPointer(name, 3);
1375:   PetscAssertPointer(value, 6);
1376:   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
1377:   PetscCall(PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value));
1378:   PetscFunctionReturn(PETSC_SUCCESS);
1379: }

1381: static PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
1382: {
1383:   htri_t exists;
1384:   hid_t  group;

1386:   PetscFunctionBegin;
1387:   PetscCallHDF5Return(exists, H5Lexists, (h5, name, H5P_DEFAULT));
1388:   if (exists) PetscCallHDF5Return(exists, H5Oexists_by_name, (h5, name, H5P_DEFAULT));
1389:   if (!exists && createGroup) {
1390:     PetscCallHDF5Return(group, H5Gcreate2, (h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1391:     PetscCallHDF5(H5Gclose, (group));
1392:     exists = PETSC_TRUE;
1393:   }
1394:   *exists_ = (PetscBool)exists;
1395:   PetscFunctionReturn(PETSC_SUCCESS);
1396: }

1398: static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
1399: {
1400:   const char rootGroupName[] = "/";
1401:   hid_t      h5;
1402:   PetscBool  exists = PETSC_FALSE;
1403:   PetscInt   i;
1404:   int        n;
1405:   char     **hierarchy;
1406:   char       buf[PETSC_MAX_PATH_LEN] = "";

1408:   PetscFunctionBegin;
1410:   if (name) PetscAssertPointer(name, 2);
1411:   else name = rootGroupName;
1412:   if (has) {
1413:     PetscAssertPointer(has, 4);
1414:     *has = PETSC_FALSE;
1415:   }
1416:   if (otype) {
1417:     PetscAssertPointer(otype, 5);
1418:     *otype = H5O_TYPE_UNKNOWN;
1419:   }
1420:   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));

1422:   /*
1423:      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
1424:      Hence, each of them needs to be tested separately:
1425:      1) whether it's a valid link
1426:      2) whether this link resolves to an object
1427:      See H5Oexists_by_name() documentation.
1428:   */
1429:   PetscCall(PetscStrToArray(name, '/', &n, &hierarchy));
1430:   if (!n) {
1431:     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1432:     if (has) *has = PETSC_TRUE;
1433:     if (otype) *otype = H5O_TYPE_GROUP;
1434:     PetscCall(PetscStrToArrayDestroy(n, hierarchy));
1435:     PetscFunctionReturn(PETSC_SUCCESS);
1436:   }
1437:   for (i = 0; i < n; i++) {
1438:     PetscCall(PetscStrlcat(buf, "/", sizeof(buf)));
1439:     PetscCall(PetscStrlcat(buf, hierarchy[i], sizeof(buf)));
1440:     PetscCall(PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists));
1441:     if (!exists) break;
1442:   }
1443:   PetscCall(PetscStrToArrayDestroy(n, hierarchy));

1445:   /* If the object exists, get its type */
1446:   if (exists && otype) {
1447:     H5O_info_t info;

1449:     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
1450:     PetscCallHDF5(H5Oget_info_by_name, (h5, name, &info, H5P_DEFAULT));
1451:     *otype = info.type;
1452:   }
1453:   if (has) *has = exists;
1454:   PetscFunctionReturn(PETSC_SUCCESS);
1455: }

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

1460:   Collective

1462:   Input Parameters:
1463: + viewer - The `PETSCVIEWERHDF5` viewer
1464: - path   - (Optional) The path relative to the pushed group

1466:   Output Parameter:
1467: . has - Flag for group existence

1469:   Level: advanced

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

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

1477: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
1478: @*/
1479: PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has)
1480: {
1481:   H5O_type_t  type;
1482:   const char *abspath;

1484:   PetscFunctionBegin;
1486:   if (path) PetscAssertPointer(path, 2);
1487:   PetscAssertPointer(has, 3);
1488:   PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath));
1489:   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
1490:   *has = (PetscBool)(type == H5O_TYPE_GROUP);
1491:   PetscCall(PetscFree(abspath));
1492:   PetscFunctionReturn(PETSC_SUCCESS);
1493: }

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

1498:   Collective

1500:   Input Parameters:
1501: + viewer - The `PETSCVIEWERHDF5` viewer
1502: - path   - The dataset path

1504:   Output Parameter:
1505: . has - Flag whether dataset exists

1507:   Level: advanced

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

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

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

1516: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasGroup()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1517: @*/
1518: PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has)
1519: {
1520:   H5O_type_t  type;
1521:   const char *abspath;

1523:   PetscFunctionBegin;
1525:   if (path) PetscAssertPointer(path, 2);
1526:   PetscAssertPointer(has, 3);
1527:   PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath));
1528:   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
1529:   *has = (PetscBool)(type == H5O_TYPE_DATASET);
1530:   PetscCall(PetscFree(abspath));
1531:   PetscFunctionReturn(PETSC_SUCCESS);
1532: }

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

1537:   Collective

1539:   Input Parameters:
1540: + viewer - The `PETSCVIEWERHDF5` viewer
1541: - obj    - The named object

1543:   Output Parameter:
1544: . has - Flag for dataset existence

1546:   Level: advanced

1548:   Notes:
1549:   If the object is unnamed, an error occurs.

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

1553: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1554: @*/
1555: PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1556: {
1557:   size_t len;

1559:   PetscFunctionBegin;
1562:   PetscAssertPointer(has, 3);
1563:   PetscCall(PetscStrlen(obj->name, &len));
1564:   PetscCheck(len, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
1565:   PetscCall(PetscViewerHDF5HasDataset(viewer, obj->name, has));
1566:   PetscFunctionReturn(PETSC_SUCCESS);
1567: }

1569: /*@C
1570:   PetscViewerHDF5HasAttribute - Check whether an attribute exists

1572:   Collective

1574:   Input Parameters:
1575: + viewer - The `PETSCVIEWERHDF5` viewer
1576: . parent - The parent dataset/group name
1577: - name   - The attribute name

1579:   Output Parameter:
1580: . has - Flag for attribute existence

1582:   Level: advanced

1584:   Note:
1585:   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.

1587: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1588: @*/
1589: PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1590: {
1591:   const char *parentAbsPath;

1593:   PetscFunctionBegin;
1595:   if (parent) PetscAssertPointer(parent, 2);
1596:   PetscAssertPointer(name, 3);
1597:   PetscAssertPointer(has, 4);
1598:   PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
1599:   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL));
1600:   if (*has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has));
1601:   PetscCall(PetscFree(parentAbsPath));
1602:   PetscFunctionReturn(PETSC_SUCCESS);
1603: }

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

1608:   Collective

1610:   Input Parameters:
1611: + viewer - The `PETSCVIEWERHDF5` viewer
1612: . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1613: - name   - The attribute name

1615:   Output Parameter:
1616: . has - Flag for attribute existence

1618:   Level: advanced

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

1624: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1625: @*/
1626: PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1627: {
1628:   PetscFunctionBegin;
1631:   PetscAssertPointer(name, 3);
1632:   PetscAssertPointer(has, 4);
1633:   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
1634:   PetscCall(PetscViewerHDF5HasAttribute(viewer, obj->name, name, has));
1635:   PetscFunctionReturn(PETSC_SUCCESS);
1636: }

1638: static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1639: {
1640:   hid_t  h5;
1641:   htri_t hhas;

1643:   PetscFunctionBegin;
1644:   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1645:   PetscCallHDF5Return(hhas, H5Aexists_by_name, (h5, parent, name, H5P_DEFAULT));
1646:   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1647:   PetscFunctionReturn(PETSC_SUCCESS);
1648: }

1650: /*
1651:   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1652:   is attached to a communicator, in this case the attribute is a PetscViewer.
1653: */
1654: PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;

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

1659:   Collective

1661:   Input Parameter:
1662: . comm - the MPI communicator to share the `PETSCVIEWERHDF5` `PetscViewer`

1664:   Options Database Key:
1665: . -viewer_hdf5_filename <name> - name of the HDF5 file

1667:   Environmental variable:
1668: . `PETSC_VIEWER_HDF5_FILENAME` - name of the HDF5 file

1670:   Level: intermediate

1672:   Note:
1673:   Unlike almost all other PETSc routines, `PETSC_VIEWER_HDF5_()` does not return
1674:   an error code.  The HDF5 `PetscViewer` is usually used in the form
1675: $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));

1677: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerCreate()`, `PetscViewerDestroy()`
1678: @*/
1679: PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1680: {
1681:   PetscErrorCode ierr;
1682:   PetscMPIInt    mpi_ierr;
1683:   PetscBool      flg;
1684:   PetscViewer    viewer;
1685:   char           fname[PETSC_MAX_PATH_LEN];
1686:   MPI_Comm       ncomm;

1688:   PetscFunctionBegin;
1689:   ierr = PetscCommDuplicate(comm, &ncomm, NULL);
1690:   if (ierr) {
1691:     ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
1692:     PetscFunctionReturn(NULL);
1693:   }
1694:   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1695:     mpi_ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Viewer_HDF5_keyval, NULL);
1696:     if (mpi_ierr) {
1697:       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
1698:       PetscFunctionReturn(NULL);
1699:     }
1700:   }
1701:   mpi_ierr = MPI_Comm_get_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void **)&viewer, (int *)&flg);
1702:   if (mpi_ierr) {
1703:     ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
1704:     PetscFunctionReturn(NULL);
1705:   }
1706:   if (!flg) { /* PetscViewer not yet created */
1707:     ierr = PetscOptionsGetenv(ncomm, "PETSC_VIEWER_HDF5_FILENAME", fname, PETSC_MAX_PATH_LEN, &flg);
1708:     if (ierr) {
1709:       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
1710:       PetscFunctionReturn(NULL);
1711:     }
1712:     if (!flg) {
1713:       ierr = PetscStrncpy(fname, "output.h5", sizeof(fname));
1714:       if (ierr) {
1715:         ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
1716:         PetscFunctionReturn(NULL);
1717:       }
1718:     }
1719:     ierr = PetscViewerHDF5Open(ncomm, fname, FILE_MODE_WRITE, &viewer);
1720:     if (ierr) {
1721:       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
1722:       PetscFunctionReturn(NULL);
1723:     }
1724:     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
1725:     if (ierr) {
1726:       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
1727:       PetscFunctionReturn(NULL);
1728:     }
1729:     mpi_ierr = MPI_Comm_set_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void *)viewer);
1730:     if (mpi_ierr) {
1731:       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
1732:       PetscFunctionReturn(NULL);
1733:     }
1734:   }
1735:   ierr = PetscCommDestroy(&ncomm);
1736:   if (ierr) {
1737:     ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
1738:     PetscFunctionReturn(NULL);
1739:   }
1740:   PetscFunctionReturn(viewer);
1741: }