Actual source code: hdf5v.c
petsc-3.14.6 2021-03-30
1: #include <petsc/private/viewerimpl.h>
2: #include <petsc/private/viewerhdf5impl.h>
3: #include <petscviewerhdf5.h>
5: static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool*, H5O_type_t*);
6: static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool*);
8: static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char objname[], char **fullpath)
9: {
10: const char *group;
11: char buf[PETSC_MAX_PATH_LEN]="";
15: PetscViewerHDF5GetGroup(viewer, &group);
16: PetscStrcat(buf, group);
17: PetscStrcat(buf, "/");
18: PetscStrcat(buf, objname);
19: PetscStrallocpy(buf, fullpath);
20: return(0);
21: }
23: static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj)
24: {
25: PetscBool has;
26: const char *group;
30: if (!obj->name) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
31: PetscViewerHDF5HasObject(viewer, obj, &has);
32: if (!has) {
33: PetscViewerHDF5GetGroup(viewer, &group);
34: SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) %s not stored in group %s", obj->name, group);
35: }
36: return(0);
37: }
39: static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v)
40: {
41: PetscErrorCode ierr;
42: PetscBool flg = PETSC_FALSE, set;
43: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;
46: PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");
47: PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);
48: PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);
49: PetscOptionsBool("-viewer_hdf5_collective","Enable collective transfer mode","PetscViewerHDF5SetCollective",flg,&flg,&set);
50: if (set) {PetscViewerHDF5SetCollective(v,flg);}
51: PetscOptionsTail();
52: return(0);
53: }
55: static PetscErrorCode PetscViewerView_HDF5(PetscViewer v,PetscViewer viewer)
56: {
57: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;
58: PetscBool flg;
59: PetscErrorCode ierr;
62: if (hdf5->filename) {
63: PetscViewerASCIIPrintf(viewer,"Filename: %s\n",hdf5->filename);
64: }
65: PetscViewerASCIIPrintf(viewer,"Vectors with blocksize 1 saved as 2D datasets: %s\n",PetscBools[hdf5->basedimension2]);
66: PetscViewerASCIIPrintf(viewer,"Enforce single precision storage: %s\n",PetscBools[hdf5->spoutput]);
67: PetscViewerHDF5GetCollective(v,&flg);
68: PetscViewerASCIIPrintf(viewer,"MPI-IO transfer mode: %s\n",flg ? "collective" : "independent");
69: return(0);
70: }
72: static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
73: {
74: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
75: PetscErrorCode ierr;
78: PetscFree(hdf5->filename);
79: if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
80: return(0);
81: }
83: static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
84: {
85: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
86: PetscErrorCode ierr;
89: PetscStackCallHDF5(H5Pclose,(hdf5->dxpl_id));
90: PetscViewerFileClose_HDF5(viewer);
91: while (hdf5->groups) {
92: PetscViewerHDF5GroupList *tmp = hdf5->groups->next;
94: PetscFree(hdf5->groups->name);
95: PetscFree(hdf5->groups);
96: hdf5->groups = tmp;
97: }
98: PetscFree(hdf5);
99: PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);
100: PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);
101: PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);
102: PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);
103: PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);
104: return(0);
105: }
107: static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
108: {
109: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
112: hdf5->btype = type;
113: return(0);
114: }
116: static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
117: {
118: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
121: *type = hdf5->btype;
122: return(0);
123: }
125: static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
126: {
127: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
130: hdf5->basedimension2 = flg;
131: return(0);
132: }
134: /*@
135: PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
136: dimension of 2.
138: Logically Collective on PetscViewer
140: Input Parameters:
141: + viewer - the PetscViewer; if it is not hdf5 then this command is ignored
142: - flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1
144: Options Database:
145: . -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
148: Notes:
149: Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
150: of one when the dimension is lower. Others think the option is crazy.
152: Level: intermediate
154: .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
156: @*/
157: PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
158: {
163: PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));
164: return(0);
165: }
167: /*@
168: PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
169: dimension of 2.
171: Logically Collective on PetscViewer
173: Input Parameter:
174: . viewer - the PetscViewer, must be of type HDF5
176: Output Parameter:
177: . flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1
179: Notes:
180: Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
181: of one when the dimension is lower. Others think the option is crazy.
183: Level: intermediate
185: .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
187: @*/
188: PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
189: {
190: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
194: *flg = hdf5->basedimension2;
195: return(0);
196: }
198: static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
199: {
200: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
203: hdf5->spoutput = flg;
204: return(0);
205: }
207: /*@
208: PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
209: compiled with double precision PetscReal.
211: Logically Collective on PetscViewer
213: Input Parameters:
214: + viewer - the PetscViewer; if it is not hdf5 then this command is ignored
215: - flg - if PETSC_TRUE the data will be written to disk with single precision
217: Options Database:
218: . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
221: Notes:
222: Setting this option does not make any difference if PETSc is compiled with single precision
223: in the first place. It does not affect reading datasets (HDF5 handle this internally).
225: Level: intermediate
227: .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
228: PetscReal
230: @*/
231: PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
232: {
237: PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));
238: return(0);
239: }
241: /*@
242: PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
243: compiled with double precision PetscReal.
245: Logically Collective on PetscViewer
247: Input Parameter:
248: . viewer - the PetscViewer, must be of type HDF5
250: Output Parameter:
251: . flg - if PETSC_TRUE the data will be written to disk with single precision
253: Notes:
254: Setting this option does not make any difference if PETSc is compiled with single precision
255: in the first place. It does not affect reading datasets (HDF5 handle this internally).
257: Level: intermediate
259: .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
260: PetscReal
262: @*/
263: PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
264: {
265: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
269: *flg = hdf5->spoutput;
270: return(0);
271: }
273: static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg)
274: {
276: /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions
277: - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */
278: #if H5_VERSION_GE(1,10,3) && defined(H5_HAVE_PARALLEL)
279: {
280: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
281: PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT));
282: }
283: #else
284: if (flg) {
286: 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");
287: }
288: #endif
289: return(0);
290: }
292: /*@
293: PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes.
295: Logically Collective; flg must contain common value
297: Input Parameters:
298: + viewer - the PetscViewer; if it is not hdf5 then this command is ignored
299: - flg - PETSC_TRUE for collective mode; PETSC_FALSE for independent mode (default)
301: Options Database:
302: . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers
304: Notes:
305: Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better.
306: 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.
308: Developer notes:
309: In the HDF5 layer, PETSC_TRUE / PETSC_FALSE means H5Pset_dxpl_mpio() is called with H5FD_MPIO_COLLECTIVE / H5FD_MPIO_INDEPENDENT, respectively.
310: This in turn means use of MPI_File_{read,write}_all / MPI_File_{read,write} in the MPI-IO layer, respectively.
311: See HDF5 documentation and MPI-IO documentation for details.
313: Level: intermediate
315: .seealso: PetscViewerHDF5GetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open()
317: @*/
318: PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg)
319: {
325: PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg));
326: return(0);
327: }
329: static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg)
330: {
331: #if defined(H5_HAVE_PARALLEL)
332: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
333: H5FD_mpio_xfer_t mode;
334: #endif
337: #if !defined(H5_HAVE_PARALLEL)
338: *flg = PETSC_FALSE;
339: #else
340: PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode));
341: *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE;
342: #endif
343: return(0);
344: }
346: /*@
347: PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes.
349: Not Collective
351: Input Parameters:
352: . viewer - the HDF5 PetscViewer
354: Output Parameters:
355: . flg - the flag
357: Level: intermediate
359: Notes:
360: 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.
361: For more details, see PetscViewerHDF5SetCollective().
363: .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open()
365: @*/
366: PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg)
367: {
374: PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));
375: return(0);
376: }
378: static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
379: {
380: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
381: hid_t plist_id;
382: PetscErrorCode ierr;
385: if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
386: if (hdf5->filename) {PetscFree(hdf5->filename);}
387: PetscStrallocpy(name, &hdf5->filename);
388: /* Set up file access property list with parallel I/O access */
389: PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
390: #if defined(H5_HAVE_PARALLEL)
391: PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL));
392: #endif
393: /* Create or open the file collectively */
394: switch (hdf5->btype) {
395: case FILE_MODE_READ:
396: PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
397: break;
398: case FILE_MODE_APPEND:
399: PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
400: break;
401: case FILE_MODE_WRITE:
402: PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
403: break;
404: default:
405: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
406: }
407: if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
408: PetscStackCallHDF5(H5Pclose,(plist_id));
409: return(0);
410: }
412: static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
413: {
414: PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
417: *name = vhdf5->filename;
418: return(0);
419: }
421: static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
422: {
423: /*
424: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
425: PetscErrorCode ierr;
426: */
429: return(0);
430: }
432: /*MC
433: PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
436: .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
437: PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
438: PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
439: PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
441: Level: beginner
442: M*/
444: PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
445: {
446: PetscViewer_HDF5 *hdf5;
447: PetscErrorCode ierr;
450: PetscNewLog(v,&hdf5);
452: v->data = (void*) hdf5;
453: v->ops->destroy = PetscViewerDestroy_HDF5;
454: v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
455: v->ops->setup = PetscViewerSetUp_HDF5;
456: v->ops->view = PetscViewerView_HDF5;
457: v->ops->flush = NULL;
458: hdf5->btype = (PetscFileMode) -1;
459: hdf5->filename = NULL;
460: hdf5->timestep = -1;
461: hdf5->groups = NULL;
463: PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER));
465: PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);
466: PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);
467: PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);
468: PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);
469: PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);
470: PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);
471: PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);
472: PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);
473: return(0);
474: }
476: /*@C
477: PetscViewerHDF5Open - Opens a file for HDF5 input/output.
479: Collective
481: Input Parameters:
482: + comm - MPI communicator
483: . name - name of file
484: - type - type of file
485: $ FILE_MODE_WRITE - create new file for binary output
486: $ FILE_MODE_READ - open existing file for binary input
487: $ FILE_MODE_APPEND - open existing file for binary output
489: Output Parameter:
490: . hdf5v - PetscViewer for HDF5 input/output to use with the specified file
492: Options Database:
493: + -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
494: - -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
496: Level: beginner
498: Note:
499: This PetscViewer should be destroyed with PetscViewerDestroy().
502: .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
503: PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
504: MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
505: @*/
506: PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
507: {
511: PetscViewerCreate(comm, hdf5v);
512: PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);
513: PetscViewerFileSetMode(*hdf5v, type);
514: PetscViewerFileSetName(*hdf5v, name);
515: PetscViewerSetFromOptions(*hdf5v);
516: return(0);
517: }
519: /*@C
520: PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
522: Not collective
524: Input Parameter:
525: . viewer - the PetscViewer
527: Output Parameter:
528: . file_id - The file id
530: Level: intermediate
532: .seealso: PetscViewerHDF5Open()
533: @*/
534: PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
535: {
536: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
540: if (file_id) *file_id = hdf5->file_id;
541: return(0);
542: }
544: /*@C
545: PetscViewerHDF5PushGroup - Set the current HDF5 group for output
547: Not collective
549: Input Parameters:
550: + viewer - the PetscViewer
551: - name - The group name
553: Level: intermediate
555: Note: The group name being NULL, empty string, or a sequence of all slashes (e.g. "///") is always internally stored as NULL and interpreted as "/".
557: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
558: @*/
559: PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[])
560: {
561: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
562: PetscViewerHDF5GroupList *groupNode;
563: PetscErrorCode ierr;
568: if (name && name[0]) {
569: size_t i,len;
570: PetscStrlen(name, &len);
571: for (i=0; i<len; i++) if (name[i] != '/') break;
572: if (i == len) name = NULL;
573: } else name = NULL;
574: PetscNew(&groupNode);
575: PetscStrallocpy(name, (char**) &groupNode->name);
576: groupNode->next = hdf5->groups;
577: hdf5->groups = groupNode;
578: return(0);
579: }
581: /*@
582: PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
584: Not collective
586: Input Parameter:
587: . viewer - the PetscViewer
589: Level: intermediate
591: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
592: @*/
593: PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer)
594: {
595: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
596: PetscViewerHDF5GroupList *groupNode;
597: PetscErrorCode ierr;
601: if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
602: groupNode = hdf5->groups;
603: hdf5->groups = hdf5->groups->next;
604: PetscFree(groupNode->name);
605: PetscFree(groupNode);
606: return(0);
607: }
609: /*@C
610: PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup().
611: If none has been assigned, returns NULL.
613: Not collective
615: Input Parameter:
616: . viewer - the PetscViewer
618: Output Parameter:
619: . name - The group name
621: Level: intermediate
623: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
624: @*/
625: PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[])
626: {
627: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
632: if (hdf5->groups) *name = hdf5->groups->name;
633: else *name = NULL;
634: return(0);
635: }
637: /*@
638: PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(),
639: and return this group's ID and file ID.
640: If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID.
642: Not collective
644: Input Parameter:
645: . viewer - the PetscViewer
647: Output Parameter:
648: + fileId - The HDF5 file ID
649: - groupId - The HDF5 group ID
651: Notes:
652: If the viewer is writable, the group is created if it doesn't exist yet.
654: Level: intermediate
656: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
657: @*/
658: PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
659: {
660: hid_t file_id;
661: H5O_type_t type;
662: const char *groupName = NULL;
663: PetscBool create;
667: PetscViewerWritable(viewer, &create);
668: PetscViewerHDF5GetFileId(viewer, &file_id);
669: PetscViewerHDF5GetGroup(viewer, &groupName);
670: PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);
671: if (type != H5O_TYPE_GROUP) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s resolves to something which is not a group", groupName);
672: PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT));
673: *fileId = file_id;
674: return(0);
675: }
677: /*@
678: PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
680: Not collective
682: Input Parameter:
683: . viewer - the PetscViewer
685: Level: intermediate
687: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
688: @*/
689: PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
690: {
691: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
695: ++hdf5->timestep;
696: return(0);
697: }
699: /*@
700: PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
701: of -1 disables blocking with timesteps.
703: Not collective
705: Input Parameters:
706: + viewer - the PetscViewer
707: - timestep - The timestep number
709: Level: intermediate
711: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
712: @*/
713: PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
714: {
715: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
719: hdf5->timestep = timestep;
720: return(0);
721: }
723: /*@
724: PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
726: Not collective
728: Input Parameter:
729: . viewer - the PetscViewer
731: Output Parameter:
732: . timestep - The timestep number
734: Level: intermediate
736: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
737: @*/
738: PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
739: {
740: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
745: *timestep = hdf5->timestep;
746: return(0);
747: }
749: /*@C
750: PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
752: Not collective
754: Input Parameter:
755: . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
757: Output Parameter:
758: . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
760: Level: advanced
762: .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
763: @*/
764: PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
765: {
767: if (ptype == PETSC_INT)
768: #if defined(PETSC_USE_64BIT_INDICES)
769: *htype = H5T_NATIVE_LLONG;
770: #else
771: *htype = H5T_NATIVE_INT;
772: #endif
773: else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE;
774: else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG;
775: else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT;
776: else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT;
777: else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT;
778: else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT;
779: else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR;
780: else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
781: else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1);
782: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
783: return(0);
784: }
786: /*@C
787: PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
789: Not collective
791: Input Parameter:
792: . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
794: Output Parameter:
795: . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
797: Level: advanced
799: .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
800: @*/
801: PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
802: {
804: #if defined(PETSC_USE_64BIT_INDICES)
805: if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG;
806: else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT;
807: #else
808: if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT;
809: #endif
810: else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
811: else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG;
812: else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT;
813: else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT;
814: else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR;
815: else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR;
816: else if (htype == H5T_C_S1) *ptype = PETSC_STRING;
817: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
818: return(0);
819: }
821: /*@C
822: PetscViewerHDF5WriteAttribute - Write an attribute
824: Input Parameters:
825: + viewer - The HDF5 viewer
826: . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
827: . name - The attribute name
828: . datatype - The attribute type
829: - value - The attribute value
831: Level: advanced
833: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
834: @*/
835: PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value)
836: {
837: char *parent;
838: hid_t h5, dataspace, obj, attribute, dtype;
839: PetscBool has;
847: PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);
848: PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);
849: PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);
850: PetscDataTypeToHDF5DataType(datatype, &dtype);
851: if (datatype == PETSC_STRING) {
852: size_t len;
853: PetscStrlen((const char *) value, &len);
854: PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
855: }
856: PetscViewerHDF5GetFileId(viewer, &h5);
857: PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
858: PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
859: if (has) {
860: PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
861: } else {
862: PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
863: }
864: PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
865: if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
866: PetscStackCallHDF5(H5Aclose,(attribute));
867: PetscStackCallHDF5(H5Oclose,(obj));
868: PetscStackCallHDF5(H5Sclose,(dataspace));
869: PetscFree(parent);
870: return(0);
871: }
873: /*@C
874: PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name
876: Input Parameters:
877: + viewer - The HDF5 viewer
878: . obj - The object whose name is used to lookup the parent dataset, relative to the current group.
879: . name - The attribute name
880: . datatype - The attribute type
881: - value - The attribute value
883: Notes:
884: This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
885: You might want to check first if it does using PetscViewerHDF5HasObject().
887: Level: advanced
889: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
890: @*/
891: PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
892: {
900: PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);
901: PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);
902: return(0);
903: }
905: /*@C
906: PetscViewerHDF5ReadAttribute - Read an attribute
908: Input Parameters:
909: + viewer - The HDF5 viewer
910: . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
911: . name - The attribute name
912: - datatype - The attribute type
914: Output Parameter:
915: . value - The attribute value
917: Notes: If the datatype is PETSC_STRING one must PetscFree() the obtained value when it is no longer needed.
919: Level: advanced
921: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
922: @*/
923: PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value)
924: {
925: char *parent;
926: hid_t h5, obj, attribute, atype, dtype;
927: PetscBool has;
935: PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);
936: PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);
937: if (has) {PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);}
938: if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name);
939: PetscDataTypeToHDF5DataType(datatype, &dtype);
940: PetscViewerHDF5GetFileId(viewer, &h5);
941: PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
942: PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
943: if (datatype == PETSC_STRING) {
944: size_t len;
945: PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
946: PetscStackCall("H5Tget_size",len = H5Tget_size(atype));
947: PetscMalloc((len+1) * sizeof(char), value);
948: PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
949: PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value));
950: } else {
951: PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
952: }
953: PetscStackCallHDF5(H5Aclose,(attribute));
954: /* H5Oclose can be used to close groups, datasets, or committed datatypes */
955: PetscStackCallHDF5(H5Oclose,(obj));
956: PetscFree(parent);
957: return(0);
958: }
960: /*@C
961: PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name
963: Input Parameters:
964: + viewer - The HDF5 viewer
965: . obj - The object whose name is used to lookup the parent dataset, relative to the current group.
966: . name - The attribute name
967: - datatype - The attribute type
969: Output Parameter:
970: . value - The attribute value
972: Notes:
973: This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
974: You might want to check first if it does using PetscViewerHDF5HasObject().
976: Level: advanced
978: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
979: @*/
980: PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value)
981: {
989: PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);
990: PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);
991: return(0);
992: }
994: PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
995: {
996: htri_t exists;
997: hid_t group;
1000: PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT));
1001: if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT));
1002: if (!exists && createGroup) {
1003: PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1004: PetscStackCallHDF5(H5Gclose,(group));
1005: exists = PETSC_TRUE;
1006: }
1007: *exists_ = (PetscBool) exists;
1008: return(0);
1009: }
1011: static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
1012: {
1013: const char rootGroupName[] = "/";
1014: hid_t h5;
1015: PetscBool exists=PETSC_FALSE;
1016: PetscInt i;
1017: int n;
1018: char **hierarchy;
1019: char buf[PETSC_MAX_PATH_LEN]="";
1025: else name = rootGroupName;
1026: if (has) {
1028: *has = PETSC_FALSE;
1029: }
1030: if (otype) {
1032: *otype = H5O_TYPE_UNKNOWN;
1033: }
1034: PetscViewerHDF5GetFileId(viewer, &h5);
1036: /*
1037: Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
1038: Hence, each of them needs to be tested separately:
1039: 1) whether it's a valid link
1040: 2) whether this link resolves to an object
1041: See H5Oexists_by_name() documentation.
1042: */
1043: PetscStrToArray(name,'/',&n,&hierarchy);
1044: if (!n) {
1045: /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1046: if (has) *has = PETSC_TRUE;
1047: if (otype) *otype = H5O_TYPE_GROUP;
1048: PetscStrToArrayDestroy(n,hierarchy);
1049: return(0);
1050: }
1051: for (i=0; i<n; i++) {
1052: PetscStrcat(buf,"/");
1053: PetscStrcat(buf,hierarchy[i]);
1054: PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);
1055: if (!exists) break;
1056: }
1057: PetscStrToArrayDestroy(n,hierarchy);
1059: /* If the object exists, get its type */
1060: if (exists && otype) {
1061: H5O_info_t info;
1063: /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
1064: PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT));
1065: *otype = info.type;
1066: }
1067: if (has) *has = exists;
1068: return(0);
1069: }
1071: /*@
1072: PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
1074: Input Parameters:
1075: . viewer - The HDF5 viewer
1077: Output Parameter:
1078: . has - Flag for group existence
1080: Notes:
1081: If the path exists but is not a group, this returns PETSC_FALSE as well.
1083: Level: advanced
1085: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup()
1086: @*/
1087: PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has)
1088: {
1089: H5O_type_t type;
1090: const char *name;
1096: PetscViewerHDF5GetGroup(viewer, &name);
1097: PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);
1098: *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE;
1099: return(0);
1100: }
1102: /*@
1103: PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1105: Input Parameters:
1106: + viewer - The HDF5 viewer
1107: - obj - The named object
1109: Output Parameter:
1110: . has - Flag for dataset existence; PETSC_FALSE for unnamed object
1112: Notes:
1113: If the path exists but is not a dataset, this returns PETSC_FALSE as well.
1115: Level: advanced
1117: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1118: @*/
1119: PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1120: {
1121: H5O_type_t type;
1122: char *path;
1129: *has = PETSC_FALSE;
1130: if (!obj->name) return(0);
1131: PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);
1132: PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);
1133: *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE;
1134: PetscFree(path);
1135: return(0);
1136: }
1138: /*@C
1139: PetscViewerHDF5HasAttribute - Check whether an attribute exists
1141: Input Parameters:
1142: + viewer - The HDF5 viewer
1143: . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
1144: - name - The attribute name
1146: Output Parameter:
1147: . has - Flag for attribute existence
1149: Level: advanced
1151: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1152: @*/
1153: PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has)
1154: {
1155: char *parent;
1163: PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);
1164: PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);
1165: if (*has) {PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);}
1166: PetscFree(parent);
1167: return(0);
1168: }
1170: /*@C
1171: PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name
1173: Input Parameters:
1174: + viewer - The HDF5 viewer
1175: . obj - The object whose name is used to lookup the parent dataset, relative to the current group.
1176: - name - The attribute name
1178: Output Parameter:
1179: . has - Flag for attribute existence
1181: Notes:
1182: This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1183: You might want to check first if it does using PetscViewerHDF5HasObject().
1185: Level: advanced
1187: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1188: @*/
1189: PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1190: {
1198: PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);
1199: PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);
1200: return(0);
1201: }
1203: static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1204: {
1205: hid_t h5;
1206: htri_t hhas;
1210: PetscViewerHDF5GetFileId(viewer, &h5);
1211: PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT));
1212: *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1213: return(0);
1214: }
1216: /*
1217: The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1218: is attached to a communicator, in this case the attribute is a PetscViewer.
1219: */
1220: PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1222: /*@C
1223: PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
1225: Collective
1227: Input Parameter:
1228: . comm - the MPI communicator to share the HDF5 PetscViewer
1230: Level: intermediate
1232: Options Database Keys:
1233: . -viewer_hdf5_filename <name>
1235: Environmental variables:
1236: . PETSC_VIEWER_HDF5_FILENAME
1238: Notes:
1239: Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
1240: an error code. The HDF5 PetscViewer is usually used in the form
1241: $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1243: .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1244: @*/
1245: PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1246: {
1248: PetscBool flg;
1249: PetscViewer viewer;
1250: char fname[PETSC_MAX_PATH_LEN];
1251: MPI_Comm ncomm;
1254: PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1255: if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1256: MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,NULL);
1257: if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1258: }
1259: MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1260: if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1261: if (!flg) { /* PetscViewer not yet created */
1262: PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
1263: if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1264: if (!flg) {
1265: PetscStrcpy(fname,"output.h5");
1266: if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1267: }
1268: PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
1269: if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1270: PetscObjectRegisterDestroy((PetscObject)viewer);
1271: if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1272: MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
1273: if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1274: }
1275: PetscCommDestroy(&ncomm);
1276: if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1277: PetscFunctionReturn(viewer);
1278: }