Actual source code: hdf5v.c
petsc-3.4.5 2014-06-29
1: #include <petsc-private/viewerimpl.h> /*I "petscsys.h" I*/
2: #include <petscviewerhdf5.h> /*I "petscviewerhdf5.h" I*/
4: typedef struct GroupList {
5: const char *name;
6: struct GroupList *next;
7: } GroupList;
9: typedef struct {
10: char *filename;
11: PetscFileMode btype;
12: hid_t file_id;
13: PetscInt timestep;
14: GroupList *groups;
15: } PetscViewer_HDF5;
19: static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
20: {
21: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
22: PetscErrorCode ierr;
25: PetscFree(hdf5->filename);
26: if (hdf5->file_id) H5Fclose(hdf5->file_id);
27: return(0);
28: }
32: PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
33: {
34: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
35: PetscErrorCode ierr;
38: PetscViewerFileClose_HDF5(viewer);
39: while (hdf5->groups) {
40: GroupList *tmp = hdf5->groups->next;
42: PetscFree(hdf5->groups->name);
43: PetscFree(hdf5->groups);
44: hdf5->groups = tmp;
45: }
46: PetscFree(hdf5);
47: PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);
48: PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);
49: return(0);
50: }
54: PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
55: {
56: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
60: hdf5->btype = type;
61: return(0);
62: }
66: PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
67: {
68: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
69: #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
70: MPI_Info info = MPI_INFO_NULL;
71: #endif
72: hid_t plist_id;
73: herr_t herr;
74: PetscErrorCode ierr;
77: PetscStrallocpy(name, &hdf5->filename);
78: /* Set up file access property list with parallel I/O access */
79: plist_id = H5Pcreate(H5P_FILE_ACCESS);
80: #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
81: herr = H5Pset_fapl_mpio(plist_id, PetscObjectComm((PetscObject)viewer), info);CHKERRQ(herr);
82: #endif
83: /* Create or open the file collectively */
84: switch (hdf5->btype) {
85: case FILE_MODE_READ:
86: hdf5->file_id = H5Fopen(name, H5F_ACC_RDONLY, plist_id);
87: break;
88: case FILE_MODE_APPEND:
89: hdf5->file_id = H5Fopen(name, H5F_ACC_RDWR, plist_id);
90: break;
91: case FILE_MODE_WRITE:
92: hdf5->file_id = H5Fcreate(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id);
93: break;
94: default:
95: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
96: }
97: if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
98: H5Pclose(plist_id);
99: return(0);
100: }
104: PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
105: {
106: PetscViewer_HDF5 *hdf5;
107: PetscErrorCode ierr;
110: PetscNewLog(v, PetscViewer_HDF5, &hdf5);
112: v->data = (void*) hdf5;
113: v->ops->destroy = PetscViewerDestroy_HDF5;
114: v->ops->flush = 0;
115: hdf5->btype = (PetscFileMode) -1;
116: hdf5->filename = 0;
117: hdf5->timestep = -1;
118: hdf5->groups = NULL;
120: PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);
121: PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);
122: return(0);
123: }
127: /*@C
128: PetscViewerHDF5Open - Opens a file for HDF5 input/output.
130: Collective on MPI_Comm
132: Input Parameters:
133: + comm - MPI communicator
134: . name - name of file
135: - type - type of file
136: $ FILE_MODE_WRITE - create new file for binary output
137: $ FILE_MODE_READ - open existing file for binary input
138: $ FILE_MODE_APPEND - open existing file for binary output
140: Output Parameter:
141: . hdf5v - PetscViewer for HDF5 input/output to use with the specified file
143: Level: beginner
145: Note:
146: This PetscViewer should be destroyed with PetscViewerDestroy().
148: Concepts: HDF5 files
149: Concepts: PetscViewerHDF5^creating
151: .seealso: PetscViewerASCIIOpen(), PetscViewerSetFormat(), PetscViewerDestroy(),
152: VecView(), MatView(), VecLoad(), MatLoad(),
153: PetscFileMode, PetscViewer
154: @*/
155: PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
156: {
160: PetscViewerCreate(comm, hdf5v);
161: PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);
162: PetscViewerFileSetMode(*hdf5v, type);
163: PetscViewerFileSetName(*hdf5v, name);
164: return(0);
165: }
169: /*@C
170: PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
172: Not collective
174: Input Parameter:
175: . viewer - the PetscViewer
177: Output Parameter:
178: . file_id - The file id
180: Level: intermediate
182: .seealso: PetscViewerHDF5Open()
183: @*/
184: PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
185: {
186: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
190: if (file_id) *file_id = hdf5->file_id;
191: return(0);
192: }
196: /*@C
197: PetscViewerHDF5PushGroup - Set the current HDF5 group for output
199: Not collective
201: Input Parameters:
202: + viewer - the PetscViewer
203: - name - The group name
205: Level: intermediate
207: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
208: @*/
209: PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
210: {
211: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
212: GroupList *groupNode;
213: PetscErrorCode ierr;
218: PetscMalloc(sizeof(GroupList), &groupNode);
219: PetscStrallocpy(name, (char**) &groupNode->name);
221: groupNode->next = hdf5->groups;
222: hdf5->groups = groupNode;
223: return(0);
224: }
228: /*@
229: PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
231: Not collective
233: Input Parameter:
234: . viewer - the PetscViewer
236: Level: intermediate
238: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup()
239: @*/
240: PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer)
241: {
242: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
243: GroupList *groupNode;
244: PetscErrorCode ierr;
248: if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
249: groupNode = hdf5->groups;
250: hdf5->groups = hdf5->groups->next;
251: PetscFree(groupNode->name);
252: PetscFree(groupNode);
253: return(0);
254: }
258: /*@C
259: PetscViewerHDF5GetGroup - Get the current HDF5 group for output. If none has been assigned, returns NULL.
261: Not collective
263: Input Parameter:
264: . viewer - the PetscViewer
266: Output Parameter:
267: . name - The group name
269: Level: intermediate
271: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup()
272: @*/
273: PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
274: {
275: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
280: if (hdf5->groups) *name = hdf5->groups->name;
281: else *name = NULL;
282: return(0);
283: }
287: /*@
288: PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
290: Not collective
292: Input Parameter:
293: . viewer - the PetscViewer
295: Level: intermediate
297: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
298: @*/
299: PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
300: {
301: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
305: ++hdf5->timestep;
306: return(0);
307: }
311: /*@
312: PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
313: of -1 disables blocking with timesteps.
315: Not collective
317: Input Parameters:
318: + viewer - the PetscViewer
319: - timestep - The timestep number
321: Level: intermediate
323: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
324: @*/
325: PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
326: {
327: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
331: hdf5->timestep = timestep;
332: return(0);
333: }
337: /*@
338: PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
340: Not collective
342: Input Parameter:
343: . viewer - the PetscViewer
345: Output Parameter:
346: . timestep - The timestep number
348: Level: intermediate
350: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
351: @*/
352: PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
353: {
354: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
359: *timestep = hdf5->timestep;
360: return(0);
361: }
363: #if defined(oldhdf4stuff)
366: PetscErrorCode PetscViewerHDF5WriteSDS(PetscViewer viewer, float *xf, int d, int *dims,int bs)
367: {
368: int i;
369: PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
370: int32 sds_id,zero32[3],dims32[3];
373: for (i = 0; i < d; i++) {
374: zero32[i] = 0;
375: dims32[i] = dims[i];
376: }
377: sds_id = SDcreate(vhdf5->sd_id, "Vec", DFNT_FLOAT32, d, dims32);
378: if (sds_id < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SDcreate failed");
379: SDwritedata(sds_id, zero32, 0, dims32, xf);
380: SDendaccess(sds_id);
381: return(0);
382: }
384: #endif