xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 2205254efee3a00a594e5e2a3a70f74dcb40bc03)
1 #include <petsc-private/viewerimpl.h>    /*I   "petscsys.h"   I*/
2 #include <hdf5.h>
3 
4 typedef struct GroupList {
5   const char       *name;
6   struct GroupList *next;
7 } GroupList;
8 
9 typedef struct {
10   char          *filename;
11   PetscFileMode btype;
12   hid_t         file_id;
13   PetscInt      timestep;
14   GroupList     *groups;
15 } PetscViewer_HDF5;
16 
17 #undef __FUNCT__
18 #define __FUNCT__ "PetscViewerFileClose_HDF5"
19 static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
20 {
21   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
22   PetscErrorCode   ierr;
23 
24   PetscFunctionBegin;
25   ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);
26   if (hdf5->file_id) H5Fclose(hdf5->file_id);
27   PetscFunctionReturn(0);
28 }
29 
30 #undef __FUNCT__
31 #define __FUNCT__ "PetscViewerDestroy_HDF5"
32 PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
33 {
34   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
35   PetscErrorCode   ierr;
36 
37   PetscFunctionBegin;
38   ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr);
39   if (hdf5->groups)
40     while (hdf5->groups) {
41       GroupList *tmp = hdf5->groups->next;
42 
43       ierr         = PetscFree(hdf5->groups->name);CHKERRQ(ierr);
44       ierr         = PetscFree(hdf5->groups);CHKERRQ(ierr);
45       hdf5->groups = tmp;
46     }
47   ierr = PetscFree(hdf5);CHKERRQ(ierr);
48   PetscFunctionReturn(0);
49 }
50 
51 EXTERN_C_BEGIN
52 #undef __FUNCT__
53 #define __FUNCT__ "PetscViewerFileSetMode_HDF5"
54 PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
55 {
56   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
57 
58   PetscFunctionBegin;
59   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
60   hdf5->btype = type;
61   PetscFunctionReturn(0);
62 }
63 EXTERN_C_END
64 
65 EXTERN_C_BEGIN
66 #undef __FUNCT__
67 #define __FUNCT__ "PetscViewerFileSetName_HDF5"
68 PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
69 {
70   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
71 #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
72   MPI_Info          info = MPI_INFO_NULL;
73 #endif
74   hid_t             plist_id;
75   herr_t            herr;
76   PetscErrorCode    ierr;
77 
78   PetscFunctionBegin;
79   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
80   /* Set up file access property list with parallel I/O access */
81   plist_id = H5Pcreate(H5P_FILE_ACCESS);
82 #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
83   herr = H5Pset_fapl_mpio(plist_id, ((PetscObject)viewer)->comm, info);CHKERRQ(herr);
84 #endif
85   /* Create or open the file collectively */
86   switch (hdf5->btype) {
87   case FILE_MODE_READ:
88     hdf5->file_id = H5Fopen(name, H5F_ACC_RDONLY, plist_id);
89     break;
90   case FILE_MODE_APPEND:
91     hdf5->file_id = H5Fopen(name, H5F_ACC_RDWR, plist_id);
92     break;
93   case FILE_MODE_WRITE:
94     hdf5->file_id = H5Fcreate(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id);
95     break;
96   default:
97     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
98   }
99   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
100   viewer->format = PETSC_VIEWER_NOFORMAT;
101   H5Pclose(plist_id);
102   PetscFunctionReturn(0);
103 }
104 EXTERN_C_END
105 
106 EXTERN_C_BEGIN
107 #undef __FUNCT__
108 #define __FUNCT__ "PetscViewerCreate_HDF5"
109 PetscErrorCode  PetscViewerCreate_HDF5(PetscViewer v)
110 {
111   PetscViewer_HDF5 *hdf5;
112   PetscErrorCode   ierr;
113 
114   PetscFunctionBegin;
115   ierr = PetscNewLog(v, PetscViewer_HDF5, &hdf5);CHKERRQ(ierr);
116 
117   v->data         = (void*) hdf5;
118   v->ops->destroy = PetscViewerDestroy_HDF5;
119   v->ops->flush   = 0;
120   v->iformat      = 0;
121   hdf5->btype     = (PetscFileMode) -1;
122   hdf5->filename  = 0;
123   hdf5->timestep  = -1;
124   hdf5->groups    = PETSC_NULL;
125 
126   ierr = PetscObjectComposeFunctionDynamic((PetscObject)v,"PetscViewerFileSetName_C","PetscViewerFileSetName_HDF5",
127                                            PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
128   ierr = PetscObjectComposeFunctionDynamic((PetscObject)v,"PetscViewerFileSetMode_C","PetscViewerFileSetMode_HDF5",
129                                            PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
130   PetscFunctionReturn(0);
131 }
132 EXTERN_C_END
133 
134 #undef __FUNCT__
135 #define __FUNCT__ "PetscViewerHDF5Open"
136 /*@C
137    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
138 
139    Collective on MPI_Comm
140 
141    Input Parameters:
142 +  comm - MPI communicator
143 .  name - name of file
144 -  type - type of file
145 $    FILE_MODE_WRITE - create new file for binary output
146 $    FILE_MODE_READ - open existing file for binary input
147 $    FILE_MODE_APPEND - open existing file for binary output
148 
149    Output Parameter:
150 .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
151 
152    Level: beginner
153 
154    Note:
155    This PetscViewer should be destroyed with PetscViewerDestroy().
156 
157    Concepts: HDF5 files
158    Concepts: PetscViewerHDF5^creating
159 
160 .seealso: PetscViewerASCIIOpen(), PetscViewerSetFormat(), PetscViewerDestroy(),
161           VecView(), MatView(), VecLoad(), MatLoad(),
162           PetscFileMode, PetscViewer
163 @*/
164 PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
165 {
166   PetscErrorCode ierr;
167 
168   PetscFunctionBegin;
169   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
170   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
171   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
172   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
173   PetscFunctionReturn(0);
174 }
175 
176 #undef __FUNCT__
177 #define __FUNCT__ "PetscViewerHDF5GetFileId"
178 /*@C
179   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
180 
181   Not collective
182 
183   Input Parameter:
184 . viewer - the PetscViewer
185 
186   Output Parameter:
187 . file_id - The file id
188 
189   Level: intermediate
190 
191 .seealso: PetscViewerHDF5Open()
192 @*/
193 PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
194 {
195   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
196 
197   PetscFunctionBegin;
198   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
199   if (file_id) *file_id = hdf5->file_id;
200   PetscFunctionReturn(0);
201 }
202 
203 #undef __FUNCT__
204 #define __FUNCT__ "PetscViewerHDF5PushGroup"
205 /*@C
206   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
207 
208   Not collective
209 
210   Input Parameters:
211 + viewer - the PetscViewer
212 - name - The group name
213 
214   Level: intermediate
215 
216 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
217 @*/
218 PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
219 {
220   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
221   GroupList        *groupNode;
222   PetscErrorCode   ierr;
223 
224   PetscFunctionBegin;
225   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
226   PetscValidCharPointer(name,2);
227   ierr = PetscMalloc(sizeof(GroupList), &groupNode);CHKERRQ(ierr);
228   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
229 
230   groupNode->next = hdf5->groups;
231   hdf5->groups    = groupNode;
232   PetscFunctionReturn(0);
233 }
234 
235 #undef __FUNCT__
236 #define __FUNCT__ "PetscViewerHDF5PopGroup"
237 /*@C
238   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
239 
240   Not collective
241 
242   Input Parameter:
243 . viewer - the PetscViewer
244 
245   Level: intermediate
246 
247 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup()
248 @*/
249 PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
250 {
251   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
252   GroupList        *groupNode;
253   PetscErrorCode   ierr;
254 
255   PetscFunctionBegin;
256   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
257   if (!hdf5->groups) SETERRQ(((PetscObject) viewer)->comm, PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
258   groupNode    = hdf5->groups;
259   hdf5->groups = hdf5->groups->next;
260   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
261   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
262   PetscFunctionReturn(0);
263 }
264 
265 #undef __FUNCT__
266 #define __FUNCT__ "PetscViewerHDF5GetGroup"
267 /*@C
268   PetscViewerHDF5GetGroup - Get the current HDF5 group for output. If none has been assigned, returns PETSC_NULL.
269 
270   Not collective
271 
272   Input Parameter:
273 . viewer - the PetscViewer
274 
275   Output Parameter:
276 . name - The group name
277 
278   Level: intermediate
279 
280 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup()
281 @*/
282 PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
283 {
284   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
285 
286   PetscFunctionBegin;
287   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
288   PetscValidCharPointer(name,2);
289   if (hdf5->groups) *name = hdf5->groups->name;
290   else *name = PETSC_NULL;
291   PetscFunctionReturn(0);
292 }
293 
294 #undef __FUNCT__
295 #define __FUNCT__ "PetscViewerHDF5IncrementTimestep"
296 /*@C
297   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
298 
299   Not collective
300 
301   Input Parameter:
302 . viewer - the PetscViewer
303 
304   Level: intermediate
305 
306 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
307 @*/
308 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
309 {
310   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
311 
312   PetscFunctionBegin;
313   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
314   ++hdf5->timestep;
315   PetscFunctionReturn(0);
316 }
317 
318 #undef __FUNCT__
319 #define __FUNCT__ "PetscViewerHDF5SetTimestep"
320 /*@C
321   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
322   of -1 disables blocking with timesteps.
323 
324   Not collective
325 
326   Input Parameters:
327 + viewer - the PetscViewer
328 - timestep - The timestep number
329 
330   Level: intermediate
331 
332 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
333 @*/
334 PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
335 {
336   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
337 
338   PetscFunctionBegin;
339   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
340   hdf5->timestep = timestep;
341   PetscFunctionReturn(0);
342 }
343 
344 #undef __FUNCT__
345 #define __FUNCT__ "PetscViewerHDF5GetTimestep"
346 /*@C
347   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
348 
349   Not collective
350 
351   Input Parameter:
352 . viewer - the PetscViewer
353 
354   Output Parameter:
355 . timestep - The timestep number
356 
357   Level: intermediate
358 
359 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
360 @*/
361 PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
362 {
363   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
364 
365   PetscFunctionBegin;
366   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
367   PetscValidPointer(timestep,2);
368   *timestep = hdf5->timestep;
369   PetscFunctionReturn(0);
370 }
371 
372 #if defined(oldhdf4stuff)
373 #undef __FUNCT__
374 #define __FUNCT__ "PetscViewerHDF5WriteSDS"
375 PetscErrorCode  PetscViewerHDF5WriteSDS(PetscViewer viewer, float *xf, int d, int *dims,int bs)
376 {
377   int              i;
378   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
379   int32            sds_id,zero32[3],dims32[3];
380 
381   PetscFunctionBegin;
382   for (i = 0; i < d; i++) {
383     zero32[i] = 0;
384     dims32[i] = dims[i];
385   }
386   sds_id = SDcreate(vhdf5->sd_id, "Vec", DFNT_FLOAT32, d, dims32);
387   if (sds_id < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SDcreate failed");
388   SDwritedata(sds_id, zero32, 0, dims32, xf);
389   SDendaccess(sds_id);
390   PetscFunctionReturn(0);
391 }
392 
393 #endif
394