Actual source code: ihdf5v.c

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

  3: static PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
  4: {
  5:   htri_t exists;
  6:   hid_t  group;

  8:   PetscFunctionBegin;
  9:   PetscCallHDF5Return(exists, H5Lexists, (h5, name, H5P_DEFAULT));
 10:   if (exists) PetscCallHDF5Return(exists, H5Oexists_by_name, (h5, name, H5P_DEFAULT));
 11:   if (!exists && createGroup) {
 12:     hid_t plist_id;
 13:     PetscCallHDF5Return(plist_id, H5Pcreate, (H5P_GROUP_CREATE));
 14:     PetscCallHDF5(H5Pset_link_creation_order, (plist_id, H5P_CRT_ORDER_TRACKED | H5P_CRT_ORDER_INDEXED));
 15:     PetscCallHDF5Return(group, H5Gcreate2, (h5, name, H5P_DEFAULT, plist_id, H5P_DEFAULT));
 16:     PetscCallHDF5(H5Pclose, (plist_id));
 17:     PetscCallHDF5(H5Gclose, (group));
 18:     exists = PETSC_TRUE;
 19:   }
 20:   *exists_ = (PetscBool)exists;
 21:   PetscFunctionReturn(PETSC_SUCCESS);
 22: }

 24: static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
 25: {
 26:   const char rootGroupName[] = "/";
 27:   hid_t      h5;
 28:   PetscBool  exists = PETSC_FALSE;
 29:   int        n;
 30:   char     **hierarchy;
 31:   char       buf[PETSC_MAX_PATH_LEN] = "";

 33:   PetscFunctionBegin;
 35:   if (name) PetscAssertPointer(name, 2);
 36:   else name = rootGroupName;
 37:   if (has) {
 38:     PetscAssertPointer(has, 4);
 39:     *has = PETSC_FALSE;
 40:   }
 41:   if (otype) {
 42:     PetscAssertPointer(otype, 5);
 43:     *otype = H5O_TYPE_UNKNOWN;
 44:   }
 45:   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));

 47:   /*
 48:      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
 49:      Hence, each of them needs to be tested separately:
 50:      1) whether it's a valid link
 51:      2) whether this link resolves to an object
 52:      See H5Oexists_by_name() documentation.
 53:   */
 54:   PetscCall(PetscStrToArray(name, '/', &n, &hierarchy));
 55:   if (!n) {
 56:     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
 57:     if (has) *has = PETSC_TRUE;
 58:     if (otype) *otype = H5O_TYPE_GROUP;
 59:     PetscCall(PetscStrToArrayDestroy(n, hierarchy));
 60:     PetscFunctionReturn(PETSC_SUCCESS);
 61:   }
 62:   for (PetscInt i = 0; i < n; i++) {
 63:     PetscCall(PetscStrlcat(buf, "/", sizeof(buf)));
 64:     PetscCall(PetscStrlcat(buf, hierarchy[i], sizeof(buf)));
 65:     PetscCall(PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists));
 66:     if (!exists) break;
 67:   }
 68:   PetscCall(PetscStrToArrayDestroy(n, hierarchy));

 70:   /* If the object exists, get its type */
 71:   if (exists && otype) {
 72:     H5O_info_t info;

 74:     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
 75:     PetscCallHDF5(H5Oget_info_by_name, (h5, name, &info, H5P_DEFAULT));
 76:     *otype = info.type;
 77:   }
 78:   if (has) *has = exists;
 79:   PetscFunctionReturn(PETSC_SUCCESS);
 80: }

 82: static PetscErrorCode PetscViewerHDF5GetGroup_Internal(PetscViewer viewer, const char *name[])
 83: {
 84:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

 86:   PetscFunctionBegin;
 88:   PetscAssertPointer(name, 2);
 89:   if (hdf5->groups) *name = hdf5->groups->name;
 90:   else *name = NULL;
 91:   PetscFunctionReturn(PETSC_SUCCESS);
 92: }

 94: static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
 95: {
 96:   hid_t  h5;
 97:   htri_t hhas;

 99:   PetscFunctionBegin;
100:   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
101:   PetscCallHDF5Return(hhas, H5Aexists_by_name, (h5, parent, name, H5P_DEFAULT));
102:   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
103:   PetscFunctionReturn(PETSC_SUCCESS);
104: }

106: static PetscErrorCode PetscViewerHDF5GetGroup_HDF5(PetscViewer viewer, const char path[], const char *abspath[])
107: {
108:   size_t      len;
109:   PetscBool   relative = PETSC_FALSE;
110:   const char *group;
111:   char        buf[PETSC_MAX_PATH_LEN] = "";

113:   PetscFunctionBegin;
114:   PetscCall(PetscViewerHDF5GetGroup_Internal(viewer, &group));
115:   PetscCall(PetscStrlen(path, &len));
116:   relative = (PetscBool)(!len || path[0] != '/');
117:   if (relative) {
118:     PetscCall(PetscStrncpy(buf, group, sizeof(buf)));
119:     if (!group || len) PetscCall(PetscStrlcat(buf, "/", sizeof(buf)));
120:     PetscCall(PetscStrlcat(buf, path, sizeof(buf)));
121:     PetscCall(PetscStrallocpy(buf, (char **)abspath));
122:   } else {
123:     PetscCall(PetscStrallocpy(path, (char **)abspath));
124:   }
125:   PetscFunctionReturn(PETSC_SUCCESS);
126: }

128: static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscViewer v, PetscOptionItems PetscOptionsObject)
129: {
130:   PetscBool         flg  = PETSC_FALSE, set;
131:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data;

133:   PetscFunctionBegin;
134:   PetscOptionsHeadBegin(PetscOptionsObject, "HDF5 PetscViewer Options");
135:   PetscCall(PetscOptionsBool("-viewer_hdf5_base_dimension2", "1d Vectors get 2 dimensions in HDF5", "PetscViewerHDF5SetBaseDimension2", hdf5->basedimension2, &hdf5->basedimension2, NULL));
136:   PetscCall(PetscOptionsBool("-viewer_hdf5_sp_output", "Force data to be written in single precision", "PetscViewerHDF5SetSPOutput", hdf5->spoutput, &hdf5->spoutput, NULL));
137:   PetscCall(PetscOptionsBool("-viewer_hdf5_collective", "Enable collective transfer mode", "PetscViewerHDF5SetCollective", flg, &flg, &set));
138:   if (set) PetscCall(PetscViewerHDF5SetCollective(v, flg));
139:   flg = PETSC_FALSE;
140:   PetscCall(PetscOptionsBool("-viewer_hdf5_default_timestepping", "Set default timestepping state", "PetscViewerHDF5SetDefaultTimestepping", flg, &flg, &set));
141:   if (set) PetscCall(PetscViewerHDF5SetDefaultTimestepping(v, flg));
142:   PetscCall(PetscOptionsBool("-viewer_hdf5_compress", "Enable compression", "PetscViewerHDF5SetCompress", hdf5->compress, &hdf5->compress, NULL));
143:   PetscOptionsHeadEnd();
144:   PetscFunctionReturn(PETSC_SUCCESS);
145: }

147: static PetscErrorCode PetscViewerView_HDF5(PetscViewer v, PetscViewer viewer)
148: {
149:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data;
150:   PetscBool         flg;

152:   PetscFunctionBegin;
153:   if (hdf5->filename) PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", hdf5->filename));
154:   PetscCall(PetscViewerASCIIPrintf(viewer, "Vectors with blocksize 1 saved as 2D datasets: %s\n", PetscBools[hdf5->basedimension2]));
155:   PetscCall(PetscViewerASCIIPrintf(viewer, "Enforce single precision storage: %s\n", PetscBools[hdf5->spoutput]));
156:   PetscCall(PetscViewerHDF5GetCollective(v, &flg));
157:   PetscCall(PetscViewerASCIIPrintf(viewer, "MPI-IO transfer mode: %s\n", flg ? "collective" : "independent"));
158:   PetscCall(PetscViewerASCIIPrintf(viewer, "Default timestepping: %s\n", PetscBools[hdf5->defTimestepping]));
159:   PetscCall(PetscViewerASCIIPrintf(viewer, "Compression: %s\n", PetscBools[hdf5->compress]));
160:   PetscFunctionReturn(PETSC_SUCCESS);
161: }

163: static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
164: {
165:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

167:   PetscFunctionBegin;
168:   PetscCall(PetscFree(hdf5->filename));
169:   if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id));
170:   PetscFunctionReturn(PETSC_SUCCESS);
171: }

173: static PetscErrorCode PetscViewerFlush_HDF5(PetscViewer viewer)
174: {
175:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

177:   PetscFunctionBegin;
178:   if (hdf5->file_id) PetscCallHDF5(H5Fflush, (hdf5->file_id, H5F_SCOPE_LOCAL));
179:   PetscFunctionReturn(PETSC_SUCCESS);
180: }

182: static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
183: {
184:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

186:   PetscFunctionBegin;
187:   PetscCallHDF5(H5Pclose, (hdf5->dxpl_id));
188:   PetscCall(PetscViewerFileClose_HDF5(viewer));
189:   while (hdf5->groups) {
190:     PetscViewerHDF5GroupList *tmp = hdf5->groups->next;

192:     PetscCall(PetscFree(hdf5->groups->name));
193:     PetscCall(PetscFree(hdf5->groups));
194:     hdf5->groups = tmp;
195:   }
196:   PetscCall(PetscFree(hdf5));
197:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetName_C", NULL));
198:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetName_C", NULL));
199:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetMode_C", NULL));
200:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetMode_C", NULL));
201:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetBaseDimension2_C", NULL));
202:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetBaseDimension2_C", NULL));
203:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetSPOutput_C", NULL));
204:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetCollective_C", NULL));
205:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetCollective_C", NULL));
206:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetDefaultTimestepping_C", NULL));
207:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetDefaultTimestepping_C", NULL));
208:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetCompress_C", NULL));
209:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetCompress_C", NULL));
210:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5WriteGroup_C", NULL));
211:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5PushGroup_C", NULL));
212:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5PopGroup_C", NULL));
213:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetGroup_C", NULL));
214:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5PushTimestepping_C", NULL));
215:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5PopTimestepping_C", NULL));
216:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5IsTimestepping_C", NULL));
217:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5IncrementTimestep_C", NULL));
218:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetTimestep_C", NULL));
219:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetTimestep_C", NULL));
220:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5WriteAttribute_C", NULL));
221:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5ReadAttribute_C", NULL));
222:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5HasGroup_C", NULL));
223:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5HasDataset_C", NULL));
224:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5HasAttribute_C", NULL));
225:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetSPOutput_C", NULL));
226:   PetscFunctionReturn(PETSC_SUCCESS);
227: }

229: static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
230: {
231:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

233:   PetscFunctionBegin;
234:   hdf5->btype = type;
235:   PetscFunctionReturn(PETSC_SUCCESS);
236: }

238: static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
239: {
240:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

242:   PetscFunctionBegin;
243:   *type = hdf5->btype;
244:   PetscFunctionReturn(PETSC_SUCCESS);
245: }

247: static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
248: {
249:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

251:   PetscFunctionBegin;
252:   hdf5->basedimension2 = flg;
253:   PetscFunctionReturn(PETSC_SUCCESS);
254: }

256: static PetscErrorCode PetscViewerHDF5GetBaseDimension2_HDF5(PetscViewer viewer, PetscBool *flg)
257: {
258:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

260:   PetscFunctionBegin;
261:   *flg = hdf5->basedimension2;
262:   PetscFunctionReturn(PETSC_SUCCESS);
263: }

265: static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
266: {
267:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

269:   PetscFunctionBegin;
270:   hdf5->spoutput = flg;
271:   PetscFunctionReturn(PETSC_SUCCESS);
272: }

274: static PetscErrorCode PetscViewerHDF5GetSPOutput_HDF5(PetscViewer viewer, PetscBool *flg)
275: {
276:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

278:   PetscFunctionBegin;
279:   *flg = hdf5->spoutput;
280:   PetscFunctionReturn(PETSC_SUCCESS);
281: }

283: static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg)
284: {
285:   PetscFunctionBegin;
286:   /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions
287:      - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */
288: #if H5_VERSION_GE(1, 10, 3) && defined(H5_HAVE_PARALLEL)
289:   {
290:     PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
291:     PetscCallHDF5(H5Pset_dxpl_mpio, (hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT));
292:   }
293: #else
294:   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"));
295: #endif
296:   PetscFunctionReturn(PETSC_SUCCESS);
297: }

299: static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg)
300: {
301: #if defined(H5_HAVE_PARALLEL)
302:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
303:   H5FD_mpio_xfer_t  mode;
304: #endif

306:   PetscFunctionBegin;
307: #if !defined(H5_HAVE_PARALLEL)
308:   *flg = PETSC_FALSE;
309: #else
310:   PetscCallHDF5(H5Pget_dxpl_mpio, (hdf5->dxpl_id, &mode));
311:   *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE;
312: #endif
313:   PetscFunctionReturn(PETSC_SUCCESS);
314: }

316: static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
317: {
318:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
319:   hid_t             plist_access_id;
320:   hid_t             plist_create_id;

322:   PetscFunctionBegin;
323:   if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id));
324:   PetscCall(PetscFree(hdf5->filename));
325:   PetscCall(PetscStrallocpy(name, &hdf5->filename));
326:   PetscCallHDF5Return(plist_create_id, H5Pcreate, (H5P_FILE_CREATE));
327:   PetscCallHDF5(H5Pset_link_creation_order, (plist_create_id, H5P_CRT_ORDER_TRACKED | H5P_CRT_ORDER_INDEXED));
328:   PetscCallHDF5Return(plist_access_id, H5Pcreate, (H5P_FILE_ACCESS));
329: #if defined(H5_HAVE_PARALLEL)
330:   PetscCallHDF5(H5Pset_fapl_mpio, (plist_access_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL));
331: #endif
332:   /* Create or open the file collectively */
333:   switch (hdf5->btype) {
334:   case FILE_MODE_READ:
335:     if (PetscDefined(USE_DEBUG)) {
336:       PetscMPIInt rank;
337:       PetscBool   flg;

339:       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
340:       if (rank == 0) {
341:         PetscCall(PetscTestFile(hdf5->filename, 'r', &flg));
342:         PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "File %s requested for reading does not exist", hdf5->filename);
343:       }
344:       PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer)));
345:     }
346:     PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDONLY, plist_access_id));
347:     break;
348:   case FILE_MODE_APPEND:
349:   case FILE_MODE_UPDATE: {
350:     PetscBool flg;
351:     PetscCall(PetscTestFile(hdf5->filename, 'r', &flg));
352:     if (flg) PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDWR, plist_access_id));
353:     else PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_EXCL, plist_create_id, plist_access_id));
354:     break;
355:   }
356:   case FILE_MODE_WRITE:
357:     PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_TRUNC, plist_create_id, plist_access_id));
358:     break;
359:   case FILE_MODE_UNDEFINED:
360:     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
361:   default:
362:     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[hdf5->btype]);
363:   }
364:   PetscCheck(hdf5->file_id >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
365:   PetscCallHDF5(H5Pclose, (plist_access_id));
366:   PetscCallHDF5(H5Pclose, (plist_create_id));
367:   PetscCall(PetscViewerHDF5ResetAttachedDMPlexStorageVersion(viewer));
368:   PetscFunctionReturn(PETSC_SUCCESS);
369: }

371: static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer, const char **name)
372: {
373:   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5 *)viewer->data;

375:   PetscFunctionBegin;
376:   *name = vhdf5->filename;
377:   PetscFunctionReturn(PETSC_SUCCESS);
378: }

380: static PetscErrorCode PetscViewerSetUp_HDF5(PETSC_UNUSED PetscViewer viewer)
381: {
382:   PetscFunctionBegin;
383:   PetscFunctionReturn(PETSC_SUCCESS);
384: }

386: static PetscErrorCode PetscViewerHDF5SetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool flg)
387: {
388:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

390:   PetscFunctionBegin;
391:   hdf5->defTimestepping = flg;
392:   PetscFunctionReturn(PETSC_SUCCESS);
393: }

395: static PetscErrorCode PetscViewerHDF5GetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool *flg)
396: {
397:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

399:   PetscFunctionBegin;
400:   *flg = hdf5->defTimestepping;
401:   PetscFunctionReturn(PETSC_SUCCESS);
402: }

404: static PetscErrorCode PetscViewerHDF5SetCompress_HDF5(PetscViewer viewer, PetscBool flg)
405: {
406:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

408:   PetscFunctionBegin;
409:   hdf5->compress = flg;
410:   PetscFunctionReturn(PETSC_SUCCESS);
411: }

413: static PetscErrorCode PetscViewerHDF5GetCompress_HDF5(PetscViewer viewer, PetscBool *flg)
414: {
415:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

417:   PetscFunctionBegin;
418:   *flg = hdf5->compress;
419:   PetscFunctionReturn(PETSC_SUCCESS);
420: }

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

425:   Not Collective

427:   Input Parameter:
428: . viewer - the `PetscViewer`

430:   Output Parameter:
431: . file_id - The file id

433:   Level: intermediate

435: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`
436: @*/
437: PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
438: {
439:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

441:   PetscFunctionBegin;
443:   if (file_id) *file_id = hdf5->file_id;
444:   PetscFunctionReturn(PETSC_SUCCESS);
445: }

447: static PetscErrorCode PetscViewerHDF5PushGroup_HDF5(PetscViewer viewer, const char name[])
448: {
449:   PetscViewer_HDF5         *hdf5 = (PetscViewer_HDF5 *)viewer->data;
450:   PetscViewerHDF5GroupList *groupNode;
451:   size_t                    i, len;
452:   char                      buf[PETSC_MAX_PATH_LEN];
453:   const char               *gname;

455:   PetscFunctionBegin;
456:   PetscCall(PetscStrlen(name, &len));
457:   gname = NULL;
458:   if (len) {
459:     if (len == 1 && name[0] == '.') {
460:       /* use current name */
461:       gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL;
462:     } else if (name[0] == '/') {
463:       /* absolute */
464:       for (i = 1; i < len; i++) {
465:         if (name[i] != '/') {
466:           gname = name;
467:           break;
468:         }
469:       }
470:     } else {
471:       /* relative */
472:       const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : "";
473:       PetscCall(PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name));
474:       gname = buf;
475:     }
476:   }
477:   PetscCall(PetscNew(&groupNode));
478:   PetscCall(PetscStrallocpy(gname, (char **)&groupNode->name));
479:   groupNode->next = hdf5->groups;
480:   hdf5->groups    = groupNode;
481:   PetscFunctionReturn(PETSC_SUCCESS);
482: }

484: static PetscErrorCode PetscViewerHDF5PopGroup_HDF5(PetscViewer viewer)
485: {
486:   PetscViewer_HDF5         *hdf5 = (PetscViewer_HDF5 *)viewer->data;
487:   PetscViewerHDF5GroupList *groupNode;

489:   PetscFunctionBegin;
490:   PetscCheck(hdf5->groups, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
491:   groupNode    = hdf5->groups;
492:   hdf5->groups = hdf5->groups->next;
493:   PetscCall(PetscFree(groupNode->name));
494:   PetscCall(PetscFree(groupNode));
495:   PetscFunctionReturn(PETSC_SUCCESS);
496: }

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

503:   Not Collective

505:   Input Parameters:
506: + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
507: - path   - (Optional) The path relative to the pushed group

509:   Output Parameters:
510: + fileId  - The HDF5 file ID
511: - groupId - The HDF5 group ID

513:   Level: intermediate

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

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

521: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5WriteGroup()`
522: @*/
523: PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, const char path[], hid_t *fileId, hid_t *groupId)
524: {
525:   hid_t       file_id;
526:   H5O_type_t  type;
527:   const char *fileName  = NULL;
528:   const char *groupName = NULL;
529:   PetscBool   writable, has;

531:   PetscFunctionBegin;
533:   if (path) PetscAssertPointer(path, 2);
534:   PetscAssertPointer(fileId, 3);
535:   PetscAssertPointer(groupId, 4);
536:   PetscCall(PetscViewerWritable(viewer, &writable));
537:   PetscCall(PetscViewerHDF5GetFileId(viewer, &file_id));
538:   PetscCall(PetscViewerFileGetName(viewer, &fileName));
539:   PetscCall(PetscViewerHDF5GetGroup(viewer, path, &groupName));
540:   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type));
541:   if (!has) {
542:     PetscCheck(writable, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Group %s does not exist and file %s is not open for writing", groupName, fileName);
543:     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName);
544:   }
545:   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);
546:   PetscCallHDF5Return(*groupId, H5Gopen2, (file_id, groupName, H5P_DEFAULT));
547:   PetscCall(PetscFree(groupName));
548:   *fileId = file_id;
549:   PetscFunctionReturn(PETSC_SUCCESS);
550: }

552: static PetscErrorCode PetscViewerHDF5WriteGroup_HDF5(PetscViewer viewer, const char path[])
553: {
554:   hid_t fileId, groupId;

556:   PetscFunctionBegin;
557:   PetscCall(PetscViewerHDF5OpenGroup(viewer, path, &fileId, &groupId)); // make sure group is actually created
558:   PetscCallHDF5(H5Gclose, (groupId));
559:   PetscFunctionReturn(PETSC_SUCCESS);
560: }

562: static PetscErrorCode PetscViewerHDF5PushTimestepping_HDF5(PetscViewer viewer)
563: {
564:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

566:   PetscFunctionBegin;
567:   PetscCheck(!hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed");
568:   hdf5->timestepping = PETSC_TRUE;
569:   if (hdf5->timestep < 0) hdf5->timestep = 0;
570:   PetscFunctionReturn(PETSC_SUCCESS);
571: }

573: static PetscErrorCode PetscViewerHDF5PopTimestepping_HDF5(PetscViewer viewer)
574: {
575:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

577:   PetscFunctionBegin;
578:   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
579:   hdf5->timestepping = PETSC_FALSE;
580:   PetscFunctionReturn(PETSC_SUCCESS);
581: }

583: static PetscErrorCode PetscViewerHDF5IsTimestepping_HDF5(PetscViewer viewer, PetscBool *flg)
584: {
585:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

587:   PetscFunctionBegin;
588:   *flg = hdf5->timestepping;
589:   PetscFunctionReturn(PETSC_SUCCESS);
590: }

592: static PetscErrorCode PetscViewerHDF5IncrementTimestep_HDF5(PetscViewer viewer)
593: {
594:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

596:   PetscFunctionBegin;
597:   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
598:   ++hdf5->timestep;
599:   PetscFunctionReturn(PETSC_SUCCESS);
600: }

602: static PetscErrorCode PetscViewerHDF5SetTimestep_HDF5(PetscViewer viewer, PetscInt timestep)
603: {
604:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

606:   PetscFunctionBegin;
607:   PetscCheck(timestep >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %" PetscInt_FMT " is negative", timestep);
608:   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
609:   hdf5->timestep = timestep;
610:   PetscFunctionReturn(PETSC_SUCCESS);
611: }

613: static PetscErrorCode PetscViewerHDF5GetTimestep_HDF5(PetscViewer viewer, PetscInt *timestep)
614: {
615:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;

617:   PetscFunctionBegin;
618:   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
619:   *timestep = hdf5->timestep;
620:   PetscFunctionReturn(PETSC_SUCCESS);
621: }

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

626:   Not Collective

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

631:   Output Parameter:
632: . htype - the HDF5 datatype

634:   Level: advanced

636: .seealso: [](sec_viewers), `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()`
637: @*/
638: PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
639: {
640:   PetscFunctionBegin;
641:   if (ptype == PETSC_INT) *htype = PetscDefined(USE_64BIT_INDICES) ? H5T_NATIVE_LLONG : H5T_NATIVE_INT;
642:   else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE;
643:   else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG;
644:   else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT;
645:   else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT;
646:   else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_HBOOL;
647:   else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT;
648:   else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR;
649:   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
650:   else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1);
651:   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
652:   PetscFunctionReturn(PETSC_SUCCESS);
653: }

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

658:   Not Collective

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

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

666:   Level: advanced

668: .seealso: [](sec_viewers), `PetscDataType`
669: @*/
670: PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
671: {
672:   PetscFunctionBegin;
673:   if (htype == H5T_NATIVE_INT) *ptype = PetscDefined(USE_64BIT_INDICES) ? PETSC_LONG : PETSC_INT;
674: #if defined(PETSC_USE_64BIT_INDICES)
675:   else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT;
676: #endif
677:   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
678:   else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG;
679:   else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT;
680:   else if (htype == H5T_NATIVE_HBOOL) *ptype = PETSC_BOOL;
681:   else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT;
682:   else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR;
683:   else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR;
684:   else if (htype == H5T_C_S1) *ptype = PETSC_STRING;
685:   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
686:   PetscFunctionReturn(PETSC_SUCCESS);
687: }

689: static PetscErrorCode PetscViewerHDF5WriteAttribute_HDF5(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
690: {
691:   const char *parentAbsPath;
692:   hid_t       h5, dataspace, obj, attribute, dtype;
693:   PetscBool   has;

695:   PetscFunctionBegin;
696:   PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
697:   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL));
698:   PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
699:   PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
700:   if (datatype == PETSC_STRING) {
701:     size_t len;
702:     PetscCall(PetscStrlen((const char *)value, &len));
703:     PetscCallHDF5(H5Tset_size, (dtype, len + 1));
704:   }
705:   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
706:   PetscCallHDF5Return(dataspace, H5Screate, (H5S_SCALAR));
707:   PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT));
708:   if (has) {
709:     PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name));
710:   } else {
711:     PetscCallHDF5Return(attribute, H5Acreate2, (obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
712:   }
713:   PetscCallHDF5(H5Awrite, (attribute, dtype, value));
714:   if (datatype == PETSC_STRING) PetscCallHDF5(H5Tclose, (dtype));
715:   PetscCallHDF5(H5Aclose, (attribute));
716:   PetscCallHDF5(H5Oclose, (obj));
717:   PetscCallHDF5(H5Sclose, (dataspace));
718:   PetscCall(PetscFree(parentAbsPath));
719:   PetscFunctionReturn(PETSC_SUCCESS);
720: }

722: static PetscErrorCode PetscViewerHDF5ReadAttribute_HDF5(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value)
723: {
724:   const char *parentAbsPath;
725:   hid_t       h5, obj, attribute, dtype;
726:   PetscBool   has;

728:   PetscFunctionBegin;
729:   PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
730:   PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
731:   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL));
732:   if (has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
733:   if (!has) {
734:     if (defaultValue) {
735:       if (defaultValue != value) {
736:         if (datatype == PETSC_STRING) {
737:           PetscCall(PetscStrallocpy(*(char **)defaultValue, (char **)value));
738:         } else {
739:           size_t len;
740:           PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (dtype));
741:           PetscCall(PetscMemcpy(value, defaultValue, len));
742:         }
743:       }
744:       PetscCall(PetscFree(parentAbsPath));
745:       PetscFunctionReturn(PETSC_SUCCESS);
746:     } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name);
747:   }
748:   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
749:   PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT));
750:   PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name));
751:   if (datatype == PETSC_STRING) {
752:     size_t len;
753:     hid_t  atype;
754:     PetscCallHDF5Return(atype, H5Aget_type, (attribute));
755:     PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (atype));
756:     PetscCall(PetscMalloc((len + 1) * sizeof(char), value));
757:     PetscCallHDF5(H5Tset_size, (dtype, len + 1));
758:     PetscCallHDF5(H5Aread, (attribute, dtype, *(char **)value));
759:   } else {
760:     PetscCallHDF5(H5Aread, (attribute, dtype, value));
761:   }
762:   PetscCallHDF5(H5Aclose, (attribute));
763:   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
764:   PetscCallHDF5(H5Oclose, (obj));
765:   PetscCall(PetscFree(parentAbsPath));
766:   PetscFunctionReturn(PETSC_SUCCESS);
767: }

769: static PetscErrorCode PetscViewerHDF5HasGroup_HDF5(PetscViewer viewer, const char path[], PetscBool *has)
770: {
771:   H5O_type_t  type;
772:   const char *abspath;

774:   PetscFunctionBegin;
775:   PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath));
776:   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
777:   *has = (PetscBool)(type == H5O_TYPE_GROUP);
778:   PetscCall(PetscFree(abspath));
779:   PetscFunctionReturn(PETSC_SUCCESS);
780: }

782: static PetscErrorCode PetscViewerHDF5HasDataset_HDF5(PetscViewer viewer, const char path[], PetscBool *has)
783: {
784:   H5O_type_t  type;
785:   const char *abspath;

787:   PetscFunctionBegin;
788:   PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath));
789:   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
790:   *has = (PetscBool)(type == H5O_TYPE_DATASET);
791:   PetscCall(PetscFree(abspath));
792:   PetscFunctionReturn(PETSC_SUCCESS);
793: }

795: static PetscErrorCode PetscViewerHDF5HasAttribute_HDF5(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
796: {
797:   const char *parentAbsPath;

799:   PetscFunctionBegin;
800:   PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
801:   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL));
802:   if (*has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has));
803:   PetscCall(PetscFree(parentAbsPath));
804:   PetscFunctionReturn(PETSC_SUCCESS);
805: }

807: /*MC
808:    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file

810:   Level: beginner

812: .seealso: [](sec_viewers), `PetscViewerHDF5Open()`, `PetscViewerStringSPrintf()`, `PetscViewerSocketOpen()`, `PetscViewerDrawOpen()`, `PETSCVIEWERSOCKET`,
813:           `PetscViewerCreate()`, `PetscViewerASCIIOpen()`, `PetscViewerBinaryOpen()`, `PETSCVIEWERBINARY`, `PETSCVIEWERDRAW`, `PETSCVIEWERSTRING`,
814:           `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`,
815:           `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`
816: M*/

818: PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
819: {
820:   PetscViewer_HDF5 *hdf5;

822:   PetscFunctionBegin;
823: #if !defined(H5_HAVE_PARALLEL)
824:   {
825:     PetscMPIInt size;
826:     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)v), &size));
827:     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)");
828:   }
829: #endif

831:   PetscCall(PetscNew(&hdf5));

833:   v->data                = (void *)hdf5;
834:   v->ops->destroy        = PetscViewerDestroy_HDF5;
835:   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
836:   v->ops->setup          = PetscViewerSetUp_HDF5;
837:   v->ops->view           = PetscViewerView_HDF5;
838:   v->ops->flush          = PetscViewerFlush_HDF5;
839:   hdf5->btype            = FILE_MODE_UNDEFINED;
840:   hdf5->filename         = NULL;
841:   hdf5->timestep         = -1;
842:   hdf5->groups           = NULL;

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

846:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_HDF5));
847:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_HDF5));
848:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_HDF5));
849:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_HDF5));
850:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetBaseDimension2_C", PetscViewerHDF5SetBaseDimension2_HDF5));
851:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetBaseDimension2_C", PetscViewerHDF5GetBaseDimension2_HDF5));
852:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetSPOutput_C", PetscViewerHDF5SetSPOutput_HDF5));
853:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetCollective_C", PetscViewerHDF5SetCollective_HDF5));
854:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetCollective_C", PetscViewerHDF5GetCollective_HDF5));
855:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetDefaultTimestepping_C", PetscViewerHDF5GetDefaultTimestepping_HDF5));
856:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetDefaultTimestepping_C", PetscViewerHDF5SetDefaultTimestepping_HDF5));
857:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetCompress_C", PetscViewerHDF5GetCompress_HDF5));
858:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetCompress_C", PetscViewerHDF5SetCompress_HDF5));
859:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5WriteGroup_C", PetscViewerHDF5WriteGroup_HDF5));
860:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5PushGroup_C", PetscViewerHDF5PushGroup_HDF5));
861:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5PopGroup_C", PetscViewerHDF5PopGroup_HDF5));
862:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetGroup_C", PetscViewerHDF5GetGroup_HDF5));
863:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5PushTimestepping_C", PetscViewerHDF5PushTimestepping_HDF5));
864:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5PopTimestepping_C", PetscViewerHDF5PopTimestepping_HDF5));
865:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5IsTimestepping_C", PetscViewerHDF5IsTimestepping_HDF5));
866:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5IncrementTimestep_C", PetscViewerHDF5IncrementTimestep_HDF5));
867:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetTimestep_C", PetscViewerHDF5SetTimestep_HDF5));
868:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetTimestep_C", PetscViewerHDF5GetTimestep_HDF5));
869:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5WriteAttribute_C", PetscViewerHDF5WriteAttribute_HDF5));
870:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5ReadAttribute_C", PetscViewerHDF5ReadAttribute_HDF5));
871:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5HasGroup_C", PetscViewerHDF5HasGroup_HDF5));
872:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5HasDataset_C", PetscViewerHDF5HasDataset_HDF5));
873:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5HasAttribute_C", PetscViewerHDF5HasAttribute_HDF5));
874:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetSPOutput_C", PetscViewerHDF5GetSPOutput_HDF5));
875:   PetscFunctionReturn(PETSC_SUCCESS);
876: }