xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 00de8ff0695ff394d09a2c60082aeaab5870b6e2)
1 #include <petsc-private/viewerimpl.h>    /*I   "petscsys.h"   I*/
2 #include <petscviewerhdf5.h>    /*I   "petscviewerhdf5.h"   I*/
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, PetscObjectComm((PetscObject)viewer), 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    = NULL;
125 
126   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C","PetscViewerFileSetName_HDF5",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
127   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C","PetscViewerFileSetMode_HDF5",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
128   PetscFunctionReturn(0);
129 }
130 EXTERN_C_END
131 
132 #undef __FUNCT__
133 #define __FUNCT__ "PetscViewerHDF5Open"
134 /*@C
135    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
136 
137    Collective on MPI_Comm
138 
139    Input Parameters:
140 +  comm - MPI communicator
141 .  name - name of file
142 -  type - type of file
143 $    FILE_MODE_WRITE - create new file for binary output
144 $    FILE_MODE_READ - open existing file for binary input
145 $    FILE_MODE_APPEND - open existing file for binary output
146 
147    Output Parameter:
148 .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
149 
150    Level: beginner
151 
152    Note:
153    This PetscViewer should be destroyed with PetscViewerDestroy().
154 
155    Concepts: HDF5 files
156    Concepts: PetscViewerHDF5^creating
157 
158 .seealso: PetscViewerASCIIOpen(), PetscViewerSetFormat(), PetscViewerDestroy(),
159           VecView(), MatView(), VecLoad(), MatLoad(),
160           PetscFileMode, PetscViewer
161 @*/
162 PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
163 {
164   PetscErrorCode ierr;
165 
166   PetscFunctionBegin;
167   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
168   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
169   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
170   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
171   PetscFunctionReturn(0);
172 }
173 
174 #undef __FUNCT__
175 #define __FUNCT__ "PetscViewerHDF5GetFileId"
176 /*@C
177   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
178 
179   Not collective
180 
181   Input Parameter:
182 . viewer - the PetscViewer
183 
184   Output Parameter:
185 . file_id - The file id
186 
187   Level: intermediate
188 
189 .seealso: PetscViewerHDF5Open()
190 @*/
191 PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
192 {
193   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
194 
195   PetscFunctionBegin;
196   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
197   if (file_id) *file_id = hdf5->file_id;
198   PetscFunctionReturn(0);
199 }
200 
201 #undef __FUNCT__
202 #define __FUNCT__ "PetscViewerHDF5PushGroup"
203 /*@C
204   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
205 
206   Not collective
207 
208   Input Parameters:
209 + viewer - the PetscViewer
210 - name - The group name
211 
212   Level: intermediate
213 
214 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
215 @*/
216 PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
217 {
218   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
219   GroupList        *groupNode;
220   PetscErrorCode   ierr;
221 
222   PetscFunctionBegin;
223   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
224   PetscValidCharPointer(name,2);
225   ierr = PetscMalloc(sizeof(GroupList), &groupNode);CHKERRQ(ierr);
226   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
227 
228   groupNode->next = hdf5->groups;
229   hdf5->groups    = groupNode;
230   PetscFunctionReturn(0);
231 }
232 
233 #undef __FUNCT__
234 #define __FUNCT__ "PetscViewerHDF5PopGroup"
235 /*@C
236   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
237 
238   Not collective
239 
240   Input Parameter:
241 . viewer - the PetscViewer
242 
243   Level: intermediate
244 
245 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup()
246 @*/
247 PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
248 {
249   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
250   GroupList        *groupNode;
251   PetscErrorCode   ierr;
252 
253   PetscFunctionBegin;
254   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
255   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
256   groupNode    = hdf5->groups;
257   hdf5->groups = hdf5->groups->next;
258   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
259   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
260   PetscFunctionReturn(0);
261 }
262 
263 #undef __FUNCT__
264 #define __FUNCT__ "PetscViewerHDF5GetGroup"
265 /*@C
266   PetscViewerHDF5GetGroup - Get the current HDF5 group for output. If none has been assigned, returns NULL.
267 
268   Not collective
269 
270   Input Parameter:
271 . viewer - the PetscViewer
272 
273   Output Parameter:
274 . name - The group name
275 
276   Level: intermediate
277 
278 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup()
279 @*/
280 PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
281 {
282   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
283 
284   PetscFunctionBegin;
285   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
286   PetscValidCharPointer(name,2);
287   if (hdf5->groups) *name = hdf5->groups->name;
288   else *name = NULL;
289   PetscFunctionReturn(0);
290 }
291 
292 #undef __FUNCT__
293 #define __FUNCT__ "PetscViewerHDF5IncrementTimestep"
294 /*@C
295   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
296 
297   Not collective
298 
299   Input Parameter:
300 . viewer - the PetscViewer
301 
302   Level: intermediate
303 
304 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
305 @*/
306 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
307 {
308   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
309 
310   PetscFunctionBegin;
311   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
312   ++hdf5->timestep;
313   PetscFunctionReturn(0);
314 }
315 
316 #undef __FUNCT__
317 #define __FUNCT__ "PetscViewerHDF5SetTimestep"
318 /*@C
319   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
320   of -1 disables blocking with timesteps.
321 
322   Not collective
323 
324   Input Parameters:
325 + viewer - the PetscViewer
326 - timestep - The timestep number
327 
328   Level: intermediate
329 
330 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
331 @*/
332 PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
333 {
334   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
335 
336   PetscFunctionBegin;
337   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
338   hdf5->timestep = timestep;
339   PetscFunctionReturn(0);
340 }
341 
342 #undef __FUNCT__
343 #define __FUNCT__ "PetscViewerHDF5GetTimestep"
344 /*@C
345   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
346 
347   Not collective
348 
349   Input Parameter:
350 . viewer - the PetscViewer
351 
352   Output Parameter:
353 . timestep - The timestep number
354 
355   Level: intermediate
356 
357 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
358 @*/
359 PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
360 {
361   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
362 
363   PetscFunctionBegin;
364   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
365   PetscValidPointer(timestep,2);
366   *timestep = hdf5->timestep;
367   PetscFunctionReturn(0);
368 }
369 
370 #if defined(oldhdf4stuff)
371 #undef __FUNCT__
372 #define __FUNCT__ "PetscViewerHDF5WriteSDS"
373 PetscErrorCode  PetscViewerHDF5WriteSDS(PetscViewer viewer, float *xf, int d, int *dims,int bs)
374 {
375   int              i;
376   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
377   int32            sds_id,zero32[3],dims32[3];
378 
379   PetscFunctionBegin;
380   for (i = 0; i < d; i++) {
381     zero32[i] = 0;
382     dims32[i] = dims[i];
383   }
384   sds_id = SDcreate(vhdf5->sd_id, "Vec", DFNT_FLOAT32, d, dims32);
385   if (sds_id < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SDcreate failed");
386   SDwritedata(sds_id, zero32, 0, dims32, xf);
387   SDendaccess(sds_id);
388   PetscFunctionReturn(0);
389 }
390 
391 #endif
392