Actual source code: hdf5v.c
petsc-3.15.0 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: case FILE_MODE_UPDATE:
400: {
401: PetscBool flg;
402: PetscTestFile(hdf5->filename, 'r', &flg);
403: if (flg) PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
404: else PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id));
405: break;
406: }
407: case FILE_MODE_WRITE:
408: PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
409: break;
410: case FILE_MODE_UNDEFINED:
411: SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
412: default:
413: SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP, "Unsupported file mode %s",PetscFileModes[hdf5->btype]);
414: }
415: if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
416: PetscStackCallHDF5(H5Pclose,(plist_id));
417: return(0);
418: }
420: static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
421: {
422: PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
425: *name = vhdf5->filename;
426: return(0);
427: }
429: static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
430: {
431: /*
432: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
433: PetscErrorCode ierr;
434: */
437: return(0);
438: }
440: /*MC
441: PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
444: .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
445: PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
446: PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
447: PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
449: Level: beginner
450: M*/
452: PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
453: {
454: PetscViewer_HDF5 *hdf5;
455: PetscErrorCode ierr;
458: PetscNewLog(v,&hdf5);
460: v->data = (void*) hdf5;
461: v->ops->destroy = PetscViewerDestroy_HDF5;
462: v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
463: v->ops->setup = PetscViewerSetUp_HDF5;
464: v->ops->view = PetscViewerView_HDF5;
465: v->ops->flush = NULL;
466: hdf5->btype = FILE_MODE_UNDEFINED;
467: hdf5->filename = NULL;
468: hdf5->timestep = -1;
469: hdf5->groups = NULL;
471: PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER));
473: PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);
474: PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);
475: PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);
476: PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);
477: PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);
478: PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);
479: PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);
480: PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);
481: return(0);
482: }
484: /*@C
485: PetscViewerHDF5Open - Opens a file for HDF5 input/output.
487: Collective
489: Input Parameters:
490: + comm - MPI communicator
491: . name - name of file
492: - type - type of file
494: Output Parameter:
495: . hdf5v - PetscViewer for HDF5 input/output to use with the specified file
497: Options Database:
498: + -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
499: - -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
501: Level: beginner
503: Notes:
504: Reading is always available, regardless of the mode. Available modes are
505: + FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY]
506: . FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC]
507: . FILE_MODE_APPEND - if file exists, keep existing contents [H5Fopen() with H5F_ACC_RDWR], else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_EXCL]
508: - FILE_MODE_UPDATE - same as FILE_MODE_APPEND
510: In case of FILE_MODE_APPEND / FILE_MODE_UPDATE, any stored object (dataset, attribute) can be selectively ovewritten if the same fully qualified name (/group/path/to/object) is specified.
512: This PetscViewer should be destroyed with PetscViewerDestroy().
515: .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
516: PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
517: MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
518: @*/
519: PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
520: {
524: PetscViewerCreate(comm, hdf5v);
525: PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);
526: PetscViewerFileSetMode(*hdf5v, type);
527: PetscViewerFileSetName(*hdf5v, name);
528: PetscViewerSetFromOptions(*hdf5v);
529: return(0);
530: }
532: /*@C
533: PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
535: Not collective
537: Input Parameter:
538: . viewer - the PetscViewer
540: Output Parameter:
541: . file_id - The file id
543: Level: intermediate
545: .seealso: PetscViewerHDF5Open()
546: @*/
547: PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
548: {
549: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
553: if (file_id) *file_id = hdf5->file_id;
554: return(0);
555: }
557: /*@C
558: PetscViewerHDF5PushGroup - Set the current HDF5 group for output
560: Not collective
562: Input Parameters:
563: + viewer - the PetscViewer
564: - name - The group name
566: Level: intermediate
568: 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 "/".
570: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
571: @*/
572: PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[])
573: {
574: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
575: PetscViewerHDF5GroupList *groupNode;
576: PetscErrorCode ierr;
581: if (name && name[0]) {
582: size_t i,len;
583: PetscStrlen(name, &len);
584: for (i=0; i<len; i++) if (name[i] != '/') break;
585: if (i == len) name = NULL;
586: } else name = NULL;
587: PetscNew(&groupNode);
588: PetscStrallocpy(name, (char**) &groupNode->name);
589: groupNode->next = hdf5->groups;
590: hdf5->groups = groupNode;
591: return(0);
592: }
594: /*@
595: PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
597: Not collective
599: Input Parameter:
600: . viewer - the PetscViewer
602: Level: intermediate
604: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
605: @*/
606: PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer)
607: {
608: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
609: PetscViewerHDF5GroupList *groupNode;
610: PetscErrorCode ierr;
614: if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
615: groupNode = hdf5->groups;
616: hdf5->groups = hdf5->groups->next;
617: PetscFree(groupNode->name);
618: PetscFree(groupNode);
619: return(0);
620: }
622: /*@C
623: PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup().
624: If none has been assigned, returns NULL.
626: Not collective
628: Input Parameter:
629: . viewer - the PetscViewer
631: Output Parameter:
632: . name - The group name
634: Level: intermediate
636: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
637: @*/
638: PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[])
639: {
640: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
645: if (hdf5->groups) *name = hdf5->groups->name;
646: else *name = NULL;
647: return(0);
648: }
650: /*@
651: PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(),
652: and return this group's ID and file ID.
653: If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID.
655: Not collective
657: Input Parameter:
658: . viewer - the PetscViewer
660: Output Parameter:
661: + fileId - The HDF5 file ID
662: - groupId - The HDF5 group ID
664: Notes:
665: If the viewer is writable, the group is created if it doesn't exist yet.
667: Level: intermediate
669: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
670: @*/
671: PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
672: {
673: hid_t file_id;
674: H5O_type_t type;
675: const char *groupName = NULL;
676: PetscBool create;
680: PetscViewerWritable(viewer, &create);
681: PetscViewerHDF5GetFileId(viewer, &file_id);
682: PetscViewerHDF5GetGroup(viewer, &groupName);
683: PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);
684: if (type != H5O_TYPE_GROUP) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s resolves to something which is not a group", groupName);
685: PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT));
686: *fileId = file_id;
687: return(0);
688: }
690: /*@
691: PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
693: Not collective
695: Input Parameter:
696: . viewer - the PetscViewer
698: Level: intermediate
700: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
701: @*/
702: PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
703: {
704: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
708: ++hdf5->timestep;
709: return(0);
710: }
712: /*@
713: PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
714: of -1 disables blocking with timesteps.
716: Not collective
718: Input Parameters:
719: + viewer - the PetscViewer
720: - timestep - The timestep number
722: Level: intermediate
724: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
725: @*/
726: PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
727: {
728: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
732: hdf5->timestep = timestep;
733: return(0);
734: }
736: /*@
737: PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
739: Not collective
741: Input Parameter:
742: . viewer - the PetscViewer
744: Output Parameter:
745: . timestep - The timestep number
747: Level: intermediate
749: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
750: @*/
751: PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
752: {
753: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
758: *timestep = hdf5->timestep;
759: return(0);
760: }
762: /*@C
763: PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
765: Not collective
767: Input Parameter:
768: . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
770: Output Parameter:
771: . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
773: Level: advanced
775: .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
776: @*/
777: PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
778: {
780: if (ptype == PETSC_INT)
781: #if defined(PETSC_USE_64BIT_INDICES)
782: *htype = H5T_NATIVE_LLONG;
783: #else
784: *htype = H5T_NATIVE_INT;
785: #endif
786: else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE;
787: else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG;
788: else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT;
789: else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT;
790: else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT;
791: else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT;
792: else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR;
793: else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
794: else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1);
795: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
796: return(0);
797: }
799: /*@C
800: PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
802: Not collective
804: Input Parameter:
805: . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
807: Output Parameter:
808: . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
810: Level: advanced
812: .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
813: @*/
814: PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
815: {
817: #if defined(PETSC_USE_64BIT_INDICES)
818: if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG;
819: else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT;
820: #else
821: if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT;
822: #endif
823: else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
824: else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG;
825: else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT;
826: else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT;
827: else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR;
828: else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR;
829: else if (htype == H5T_C_S1) *ptype = PETSC_STRING;
830: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
831: return(0);
832: }
834: /*@C
835: PetscViewerHDF5WriteAttribute - Write an attribute
837: Input Parameters:
838: + viewer - The HDF5 viewer
839: . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
840: . name - The attribute name
841: . datatype - The attribute type
842: - value - The attribute value
844: Level: advanced
846: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
847: @*/
848: PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value)
849: {
850: char *parent;
851: hid_t h5, dataspace, obj, attribute, dtype;
852: PetscBool has;
860: PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);
861: PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);
862: PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);
863: PetscDataTypeToHDF5DataType(datatype, &dtype);
864: if (datatype == PETSC_STRING) {
865: size_t len;
866: PetscStrlen((const char *) value, &len);
867: PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
868: }
869: PetscViewerHDF5GetFileId(viewer, &h5);
870: PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
871: PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
872: if (has) {
873: PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
874: } else {
875: PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
876: }
877: PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
878: if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
879: PetscStackCallHDF5(H5Aclose,(attribute));
880: PetscStackCallHDF5(H5Oclose,(obj));
881: PetscStackCallHDF5(H5Sclose,(dataspace));
882: PetscFree(parent);
883: return(0);
884: }
886: /*@C
887: PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name
889: Input Parameters:
890: + viewer - The HDF5 viewer
891: . obj - The object whose name is used to lookup the parent dataset, relative to the current group.
892: . name - The attribute name
893: . datatype - The attribute type
894: - value - The attribute value
896: Notes:
897: This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
898: You might want to check first if it does using PetscViewerHDF5HasObject().
900: Level: advanced
902: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
903: @*/
904: PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
905: {
913: PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);
914: PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);
915: return(0);
916: }
918: /*@C
919: PetscViewerHDF5ReadAttribute - Read an attribute
921: Input Parameters:
922: + viewer - The HDF5 viewer
923: . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
924: . name - The attribute name
925: - datatype - The attribute type
927: Output Parameter:
928: . value - The attribute value
930: Notes: If the datatype is PETSC_STRING one must PetscFree() the obtained value when it is no longer needed.
932: Level: advanced
934: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
935: @*/
936: PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value)
937: {
938: char *parent;
939: hid_t h5, obj, attribute, atype, dtype;
940: PetscBool has;
948: PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);
949: PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);
950: if (has) {PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);}
951: if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name);
952: PetscDataTypeToHDF5DataType(datatype, &dtype);
953: PetscViewerHDF5GetFileId(viewer, &h5);
954: PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
955: PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
956: if (datatype == PETSC_STRING) {
957: size_t len;
958: PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
959: PetscStackCall("H5Tget_size",len = H5Tget_size(atype));
960: PetscMalloc((len+1) * sizeof(char), value);
961: PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
962: PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value));
963: } else {
964: PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
965: }
966: PetscStackCallHDF5(H5Aclose,(attribute));
967: /* H5Oclose can be used to close groups, datasets, or committed datatypes */
968: PetscStackCallHDF5(H5Oclose,(obj));
969: PetscFree(parent);
970: return(0);
971: }
973: /*@C
974: PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name
976: Input Parameters:
977: + viewer - The HDF5 viewer
978: . obj - The object whose name is used to lookup the parent dataset, relative to the current group.
979: . name - The attribute name
980: - datatype - The attribute type
982: Output Parameter:
983: . value - The attribute value
985: Notes:
986: This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
987: You might want to check first if it does using PetscViewerHDF5HasObject().
989: Level: advanced
991: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
992: @*/
993: PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value)
994: {
1002: PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);
1003: PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);
1004: return(0);
1005: }
1007: PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
1008: {
1009: htri_t exists;
1010: hid_t group;
1013: PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT));
1014: if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT));
1015: if (!exists && createGroup) {
1016: PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1017: PetscStackCallHDF5(H5Gclose,(group));
1018: exists = PETSC_TRUE;
1019: }
1020: *exists_ = (PetscBool) exists;
1021: return(0);
1022: }
1024: static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
1025: {
1026: const char rootGroupName[] = "/";
1027: hid_t h5;
1028: PetscBool exists=PETSC_FALSE;
1029: PetscInt i;
1030: int n;
1031: char **hierarchy;
1032: char buf[PETSC_MAX_PATH_LEN]="";
1038: else name = rootGroupName;
1039: if (has) {
1041: *has = PETSC_FALSE;
1042: }
1043: if (otype) {
1045: *otype = H5O_TYPE_UNKNOWN;
1046: }
1047: PetscViewerHDF5GetFileId(viewer, &h5);
1049: /*
1050: Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
1051: Hence, each of them needs to be tested separately:
1052: 1) whether it's a valid link
1053: 2) whether this link resolves to an object
1054: See H5Oexists_by_name() documentation.
1055: */
1056: PetscStrToArray(name,'/',&n,&hierarchy);
1057: if (!n) {
1058: /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1059: if (has) *has = PETSC_TRUE;
1060: if (otype) *otype = H5O_TYPE_GROUP;
1061: PetscStrToArrayDestroy(n,hierarchy);
1062: return(0);
1063: }
1064: for (i=0; i<n; i++) {
1065: PetscStrcat(buf,"/");
1066: PetscStrcat(buf,hierarchy[i]);
1067: PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);
1068: if (!exists) break;
1069: }
1070: PetscStrToArrayDestroy(n,hierarchy);
1072: /* If the object exists, get its type */
1073: if (exists && otype) {
1074: H5O_info_t info;
1076: /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
1077: PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT));
1078: *otype = info.type;
1079: }
1080: if (has) *has = exists;
1081: return(0);
1082: }
1084: /*@
1085: PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
1087: Input Parameters:
1088: . viewer - The HDF5 viewer
1090: Output Parameter:
1091: . has - Flag for group existence
1093: Notes:
1094: If the path exists but is not a group, this returns PETSC_FALSE as well.
1096: Level: advanced
1098: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup()
1099: @*/
1100: PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has)
1101: {
1102: H5O_type_t type;
1103: const char *name;
1109: PetscViewerHDF5GetGroup(viewer, &name);
1110: PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);
1111: *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE;
1112: return(0);
1113: }
1115: /*@
1116: PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1118: Input Parameters:
1119: + viewer - The HDF5 viewer
1120: - obj - The named object
1122: Output Parameter:
1123: . has - Flag for dataset existence; PETSC_FALSE for unnamed object
1125: Notes:
1126: If the path exists but is not a dataset, this returns PETSC_FALSE as well.
1128: Level: advanced
1130: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1131: @*/
1132: PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1133: {
1134: H5O_type_t type;
1135: char *path;
1142: *has = PETSC_FALSE;
1143: if (!obj->name) return(0);
1144: PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);
1145: PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);
1146: *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE;
1147: PetscFree(path);
1148: return(0);
1149: }
1151: /*@C
1152: PetscViewerHDF5HasAttribute - Check whether an attribute exists
1154: Input Parameters:
1155: + viewer - The HDF5 viewer
1156: . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
1157: - name - The attribute name
1159: Output Parameter:
1160: . has - Flag for attribute existence
1162: Level: advanced
1164: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1165: @*/
1166: PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has)
1167: {
1168: char *parent;
1176: PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);
1177: PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);
1178: if (*has) {PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);}
1179: PetscFree(parent);
1180: return(0);
1181: }
1183: /*@C
1184: PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name
1186: Input Parameters:
1187: + viewer - The HDF5 viewer
1188: . obj - The object whose name is used to lookup the parent dataset, relative to the current group.
1189: - name - The attribute name
1191: Output Parameter:
1192: . has - Flag for attribute existence
1194: Notes:
1195: This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1196: You might want to check first if it does using PetscViewerHDF5HasObject().
1198: Level: advanced
1200: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1201: @*/
1202: PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1203: {
1211: PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);
1212: PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);
1213: return(0);
1214: }
1216: static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1217: {
1218: hid_t h5;
1219: htri_t hhas;
1223: PetscViewerHDF5GetFileId(viewer, &h5);
1224: PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT));
1225: *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1226: return(0);
1227: }
1229: /*
1230: The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1231: is attached to a communicator, in this case the attribute is a PetscViewer.
1232: */
1233: PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1235: /*@C
1236: PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
1238: Collective
1240: Input Parameter:
1241: . comm - the MPI communicator to share the HDF5 PetscViewer
1243: Level: intermediate
1245: Options Database Keys:
1246: . -viewer_hdf5_filename <name>
1248: Environmental variables:
1249: . PETSC_VIEWER_HDF5_FILENAME
1251: Notes:
1252: Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
1253: an error code. The HDF5 PetscViewer is usually used in the form
1254: $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1256: .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1257: @*/
1258: PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1259: {
1261: PetscBool flg;
1262: PetscViewer viewer;
1263: char fname[PETSC_MAX_PATH_LEN];
1264: MPI_Comm ncomm;
1267: PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1268: if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1269: MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,NULL);
1270: if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1271: }
1272: MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1273: if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1274: if (!flg) { /* PetscViewer not yet created */
1275: PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
1276: if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1277: if (!flg) {
1278: PetscStrcpy(fname,"output.h5");
1279: if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1280: }
1281: PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
1282: if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1283: PetscObjectRegisterDestroy((PetscObject)viewer);
1284: if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1285: MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
1286: if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1287: }
1288: PetscCommDestroy(&ncomm);
1289: if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1290: PetscFunctionReturn(viewer);
1291: }