Line data Source code
1 : !--------------------------------------------------------------------------------------------------!
2 : ! CP2K: A general program to perform molecular dynamics simulations !
3 : ! Copyright 2000-2026 CP2K developers group <https://cp2k.org> !
4 : ! !
5 : ! SPDX-License-Identifier: GPL-2.0-or-later !
6 : !--------------------------------------------------------------------------------------------------!
7 :
8 : ! **************************************************************************************************
9 : !> \brief A wrapper around the HDF5 Fortran API
10 : !> \par History
11 : !> 04.2023 created [SB]
12 : !> \author Stefano Battaglia
13 : ! **************************************************************************************************
14 : MODULE hdf5_wrapper
15 :
16 : #ifdef __HDF5
17 : USE hdf5, ONLY: &
18 : h5aclose_f, h5acreate_f, h5aopen_f, h5aread_f, h5awrite_f, h5close_f, h5dclose_f, &
19 : h5dcreate_f, h5dget_space_f, h5dopen_f, h5dread_f, h5dwrite_f, h5f_acc_rdonly_f, &
20 : h5f_acc_trunc_f, h5fclose_f, h5fcreate_f, h5fopen_f, h5gclose_f, h5gcreate_f, h5gopen_f, &
21 : h5open_f, h5s_scalar_f, h5sclose_f, h5screate_f, h5screate_simple_f, &
22 : h5sget_simple_extent_npoints_f, h5t_c_s1, h5t_cset_utf8_f, h5t_enum_f, h5t_native_double, &
23 : h5t_native_integer, h5t_str_nullpad_f, h5t_string, h5tclose_f, h5tcopy_f, h5tcreate_f, &
24 : h5tenum_insert_f, h5tset_cset_f, h5tset_size_f, h5tset_strpad_f, hid_t, hsize_t, size_t
25 : #endif
26 : USE iso_c_binding, ONLY: C_LOC, &
27 : c_ptr
28 : USE cp_log_handling, ONLY: cp_logger_get_default_io_unit
29 : USE kinds, ONLY: dp
30 : #include "./base/base_uses.f90"
31 :
32 : IMPLICIT NONE
33 :
34 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'hdf5_wrapper'
35 : #ifdef __HDF5
36 : INTEGER, PARAMETER, PUBLIC :: hdf5_id = hid_t
37 : #endif
38 :
39 : CONTAINS
40 :
41 : #ifdef __HDF5
42 : ! **************************************************************************************************
43 : !> \brief Initialize the HDF5 fortran API
44 : ! **************************************************************************************************
45 0 : SUBROUTINE h5open()
46 : INTEGER :: error
47 :
48 0 : CALL h5open_f(error)
49 0 : IF (error < 0) CPABORT('ERROR: failed to initialize HDF5 interface')
50 :
51 0 : END SUBROUTINE h5open
52 :
53 : ! **************************************************************************************************
54 : !> \brief Close the HDF5 fortran API
55 : ! **************************************************************************************************
56 0 : SUBROUTINE h5close()
57 : INTEGER :: error
58 :
59 0 : CALL h5close_f(error)
60 0 : IF (error < 0) CPABORT('ERROR: failed to close HDF5 interface')
61 :
62 0 : END SUBROUTINE h5close
63 :
64 : ! **************************************************************************************************
65 : !> \brief Create a HDF5 file
66 : !> \param filename the name of the hdf5 file
67 : !> \param file_id the file id of the hdf5 file
68 : ! **************************************************************************************************
69 0 : SUBROUTINE h5fcreate(filename, file_id)
70 : CHARACTER(LEN=*), INTENT(IN) :: filename
71 : INTEGER(KIND=hid_t), INTENT(OUT) :: file_id
72 :
73 : INTEGER :: error
74 :
75 0 : CALL h5fcreate_f(filename, h5f_acc_trunc_f, file_id, error)
76 0 : IF (error < 0) CPABORT('ERROR: failed to create HDF5 file')
77 :
78 0 : END SUBROUTINE h5fcreate
79 :
80 : ! **************************************************************************************************
81 : !> \brief Open a HDF5 file
82 : !> \param filename the name of the hdf5 file
83 : !> \param file_id the file id of the hdf5 file
84 : ! **************************************************************************************************
85 0 : SUBROUTINE h5fopen(filename, file_id)
86 : CHARACTER(LEN=*), INTENT(IN) :: filename
87 : INTEGER(KIND=hid_t), INTENT(OUT) :: file_id
88 :
89 : INTEGER :: error
90 :
91 0 : CALL h5fopen_f(TRIM(filename), h5f_acc_rdonly_f, file_id, error)
92 0 : IF (error < 0) CPABORT('ERROR: failed to open HDF5 file')
93 :
94 0 : END SUBROUTINE h5fopen
95 :
96 : ! **************************************************************************************************
97 : !> \brief Close a HDF5 file
98 : !> \param file_id the file id of the hdf5 file
99 : ! **************************************************************************************************
100 0 : SUBROUTINE h5fclose(file_id)
101 : INTEGER(KIND=hid_t), INTENT(IN) :: file_id
102 :
103 : INTEGER :: error
104 :
105 0 : CALL h5fclose_f(file_id, error)
106 0 : IF (error < 0) CPABORT('ERROR: failed to close HDF5 file')
107 :
108 0 : END SUBROUTINE h5fclose
109 :
110 : ! **************************************************************************************************
111 : !> \brief Create a HDF5 group
112 : !> \param loc_id file or group identifier
113 : !> \param name name of the group
114 : !> \param grp_id group identifier
115 : ! **************************************************************************************************
116 0 : SUBROUTINE h5gcreate(loc_id, name, grp_id)
117 : INTEGER(KIND=hid_t), INTENT(IN) :: loc_id
118 : CHARACTER(LEN=*), INTENT(IN) :: name
119 : INTEGER(KIND=hid_t), INTENT(OUT) :: grp_id
120 :
121 : INTEGER :: error
122 :
123 0 : CALL h5gcreate_f(loc_id, name, grp_id, error)
124 0 : IF (error < 0) CPABORT('ERROR: failed to create HDF5 group')
125 :
126 0 : END SUBROUTINE h5gcreate
127 :
128 : ! **************************************************************************************************
129 : !> \brief Open a HDF5 group
130 : !> \param loc_id file or group identifier
131 : !> \param name name of the group
132 : !> \param grp_id group identifier
133 : ! **************************************************************************************************
134 0 : SUBROUTINE h5gopen(loc_id, name, grp_id)
135 : INTEGER(KIND=hid_t), INTENT(IN) :: loc_id
136 : CHARACTER(LEN=*), INTENT(IN) :: name
137 : INTEGER(KIND=hid_t), INTENT(OUT) :: grp_id
138 :
139 : INTEGER :: error
140 :
141 0 : CALL h5gopen_f(loc_id, name, grp_id, error)
142 0 : IF (error < 0) CPABORT('ERROR: failed to open HDF5 group')
143 :
144 0 : END SUBROUTINE h5gopen
145 :
146 : ! **************************************************************************************************
147 : !> \brief Close a HDF5 group
148 : !> \param grp_id group identifier
149 : ! **************************************************************************************************
150 0 : SUBROUTINE h5gclose(grp_id)
151 : INTEGER(KIND=hid_t), INTENT(IN) :: grp_id
152 :
153 : INTEGER :: error
154 :
155 0 : CALL h5gclose_f(grp_id, error)
156 0 : IF (error < 0) CPABORT('ERROR: failed to close HDF5 group')
157 :
158 0 : END SUBROUTINE h5gclose
159 :
160 : ! **************************************************************************************************
161 : !> \brief Write a variable-length string attribute
162 : !> \param loc_id either file id or group id
163 : !> \param attr_name the name of the attribute
164 : !> \param attr_data the attribute data, i.e. the string to write
165 : ! **************************************************************************************************
166 0 : SUBROUTINE h5awrite_varlen_string(loc_id, attr_name, attr_data)
167 : INTEGER(KIND=hid_t), INTENT(IN) :: loc_id
168 : CHARACTER(LEN=*), INTENT(IN) :: attr_name
169 : CHARACTER(LEN=*), INTENT(IN), TARGET :: attr_data
170 :
171 : INTEGER :: error, output_unit
172 : INTEGER(KIND=hid_t) :: attr_id, space_id, type_id
173 : TYPE(c_ptr) :: buffer
174 : TYPE(c_ptr), TARGET :: in_between_ptr
175 :
176 0 : output_unit = cp_logger_get_default_io_unit()
177 :
178 : ! create a scalar dataspace
179 0 : CALL h5screate_f(h5s_scalar_f, space_id, error)
180 0 : IF (error < 0) THEN
181 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
182 0 : ' ERROR: failed to create HDF5 dataspace'
183 0 : RETURN
184 : END IF
185 :
186 : ! create a variable-length string type
187 0 : CALL h5tcopy_f(h5t_string, type_id, error)
188 0 : CALL h5tset_cset_f(type_id, h5t_cset_utf8_f, error)
189 0 : CALL h5tset_strpad_f(type_id, h5t_str_nullpad_f, error)
190 :
191 : ! create the attribute
192 0 : CALL h5acreate_f(loc_id, attr_name, type_id, space_id, attr_id, error)
193 0 : IF (error < 0) THEN
194 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
195 0 : ' ERROR: failed to create HDF5 attribute'
196 0 : RETURN
197 : END IF
198 :
199 : ! weird in-between pointer needed for variable-length
200 : ! string to a scalar dataspace
201 0 : in_between_ptr = C_LOC(attr_data)
202 : ! the actual pointer to be passed
203 0 : buffer = C_LOC(in_between_ptr)
204 :
205 : ! write the string attribute to file
206 0 : CALL h5awrite_f(attr_id, type_id, buffer, error)
207 0 : IF (error < 0) THEN
208 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
209 0 : ' ERROR: failed to write HDF5 attribute'
210 0 : RETURN
211 : END IF
212 :
213 : ! close attribute
214 0 : CALL h5aclose_f(attr_id, error)
215 0 : IF (error < 0) THEN
216 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
217 0 : ' ERROR: failed to close HDF5 attribute'
218 0 : RETURN
219 : END IF
220 :
221 : ! close dataspace
222 0 : CALL h5sclose_f(space_id, error)
223 0 : IF (error < 0) THEN
224 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
225 0 : ' ERROR: failed to close HDF5 dataspace'
226 0 : RETURN
227 : END IF
228 :
229 : ! close datatype
230 0 : CALL h5tclose_f(type_id, error)
231 0 : IF (error < 0) THEN
232 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
233 0 : ' ERROR: failed to close HDF5 datatype'
234 0 : RETURN
235 : END IF
236 :
237 : END SUBROUTINE h5awrite_varlen_string
238 :
239 : ! **************************************************************************************************
240 : !> \brief Write a fixed-length string attribute
241 : !> \param loc_id either file id or group id
242 : !> \param attr_name the name of the attribute
243 : !> \param attr_data the attribute data, i.e. the string to write
244 : ! **************************************************************************************************
245 0 : SUBROUTINE h5awrite_fixlen_string(loc_id, attr_name, attr_data)
246 : INTEGER(KIND=hid_t), INTENT(IN) :: loc_id
247 : CHARACTER(LEN=*), INTENT(IN) :: attr_name
248 : CHARACTER(LEN=*), INTENT(IN), TARGET :: attr_data
249 :
250 : INTEGER :: error, output_unit
251 : INTEGER(KIND=hid_t) :: attr_id, space_id, type_id
252 : TYPE(c_ptr) :: buffer
253 :
254 0 : output_unit = cp_logger_get_default_io_unit()
255 :
256 : ! create a scalar dataspace
257 0 : CALL h5screate_f(h5s_scalar_f, space_id, error)
258 0 : IF (error < 0) THEN
259 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
260 0 : ' ERROR: failed to create HDF5 dataspace'
261 0 : RETURN
262 : END IF
263 :
264 : ! create a fixed-length string datatype
265 0 : CALL h5tcopy_f(h5t_c_s1, type_id, error)
266 0 : CALL h5tset_cset_f(type_id, h5t_cset_utf8_f, error)
267 0 : CALL h5tset_size_f(type_id, LEN(attr_data, size_t), error)
268 :
269 : ! create the attribute
270 0 : CALL h5acreate_f(loc_id, attr_name, type_id, space_id, attr_id, error)
271 0 : IF (error < 0) THEN
272 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
273 0 : ' ERROR: failed to create HDF5 attribute'
274 0 : RETURN
275 : END IF
276 :
277 : ! the actual pointer to be passed
278 0 : buffer = C_LOC(attr_data)
279 :
280 : ! write the string attribute to file
281 0 : CALL h5awrite_f(attr_id, type_id, buffer, error)
282 0 : IF (error < 0) THEN
283 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
284 0 : ' ERROR: failed to write HDF5 attribute'
285 0 : RETURN
286 : END IF
287 :
288 : ! close attribute
289 0 : CALL h5aclose_f(attr_id, error)
290 0 : IF (error < 0) THEN
291 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
292 0 : ' ERROR: failed to close HDF5 attribute'
293 0 : RETURN
294 : END IF
295 :
296 : ! close dataspace
297 0 : CALL h5sclose_f(space_id, error)
298 0 : IF (error < 0) THEN
299 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
300 0 : ' ERROR: failed to close HDF5 dataspace'
301 0 : RETURN
302 : END IF
303 :
304 : ! close datatype
305 0 : CALL h5tclose_f(type_id, error)
306 0 : IF (error < 0) THEN
307 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
308 0 : ' ERROR: failed to close HDF5 datatype'
309 0 : RETURN
310 : END IF
311 :
312 : END SUBROUTINE h5awrite_fixlen_string
313 :
314 : ! **************************************************************************************************
315 : !> \brief Write a boolean attribute
316 : !> \param loc_id either file id or group id
317 : !> \param attr_name the name of the attribute
318 : !> \param attr_data the attribute data, i.e. the logical to write (.true. or .false.)
319 : ! **************************************************************************************************
320 0 : SUBROUTINE h5awrite_boolean(loc_id, attr_name, attr_data)
321 : INTEGER(KIND=hid_t), INTENT(IN) :: loc_id
322 : CHARACTER(LEN=*), INTENT(IN) :: attr_name
323 : LOGICAL, INTENT(IN) :: attr_data
324 :
325 : INTEGER :: error, output_unit
326 : INTEGER(KIND=hid_t) :: attr_id, space_id, type_id
327 : INTEGER, TARGET :: attr_data_to_int
328 : TYPE(c_ptr) :: buffer
329 :
330 0 : output_unit = cp_logger_get_default_io_unit()
331 :
332 : ! 8-bit integers in enum bool_type
333 :
334 : ! create a scalar dataspace
335 0 : CALL h5screate_f(h5s_scalar_f, space_id, error)
336 0 : IF (error < 0) THEN
337 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
338 0 : ' ERROR: failed to create HDF5 dataspace'
339 0 : RETURN
340 : END IF
341 :
342 : ! create the datatype
343 0 : CALL h5tcreate_f(h5t_enum_f, INT(1, size_t), type_id, error)
344 0 : CALL h5tenum_insert_f(type_id, "FALSE", 0, error)
345 0 : CALL h5tenum_insert_f(type_id, "TRUE", 1, error)
346 :
347 0 : IF (attr_data) THEN
348 0 : attr_data_to_int = 1
349 : ELSE
350 0 : attr_data_to_int = 0
351 : END IF
352 : ! the C pointer to the actual data
353 0 : buffer = C_LOC(attr_data_to_int)
354 :
355 : ! create the attribute
356 0 : CALL h5acreate_f(loc_id, attr_name, type_id, space_id, attr_id, error)
357 0 : IF (error < 0) THEN
358 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
359 0 : ' ERROR: failed to create HDF5 attribute'
360 0 : RETURN
361 : END IF
362 :
363 : ! write the string attribute to file
364 0 : CALL h5awrite_f(attr_id, type_id, buffer, error)
365 0 : IF (error < 0) THEN
366 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
367 0 : ' ERROR: failed to write HDF5 attribute'
368 0 : RETURN
369 : END IF
370 :
371 : ! close attribute
372 0 : CALL h5aclose_f(attr_id, error)
373 0 : IF (error < 0) THEN
374 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
375 0 : ' ERROR: failed to close HDF5 attribute'
376 0 : RETURN
377 : END IF
378 :
379 : ! close dataspace
380 0 : CALL h5sclose_f(space_id, error)
381 0 : IF (error < 0) THEN
382 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
383 0 : ' ERROR: failed to close HDF5 dataspace'
384 0 : RETURN
385 : END IF
386 :
387 : ! close datatype
388 0 : CALL h5tclose_f(type_id, error)
389 0 : IF (error < 0) THEN
390 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
391 0 : ' ERROR: failed to close HDF5 datatype'
392 0 : RETURN
393 : END IF
394 :
395 : END SUBROUTINE h5awrite_boolean
396 :
397 : ! **************************************************************************************************
398 : !> \brief Write a (scalar) integer attribute
399 : !> \param loc_id either file id or group id
400 : !> \param attr_name the name of the attribute
401 : !> \param attr_data the attribute data, i.e. the integer to write
402 : ! **************************************************************************************************
403 0 : SUBROUTINE h5awrite_integer_scalar(loc_id, attr_name, attr_data)
404 : INTEGER(KIND=hid_t), INTENT(IN) :: loc_id
405 : CHARACTER(LEN=*), INTENT(IN) :: attr_name
406 : INTEGER, INTENT(IN), TARGET :: attr_data
407 :
408 : INTEGER :: error, output_unit
409 : INTEGER(KIND=hid_t) :: attr_id, space_id, type_id
410 : TYPE(c_ptr) :: buffer
411 :
412 0 : output_unit = cp_logger_get_default_io_unit()
413 :
414 : ! create a scalar dataspace
415 0 : CALL h5screate_f(h5s_scalar_f, space_id, error)
416 0 : IF (error < 0) THEN
417 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
418 0 : ' ERROR: failed to create HDF5 dataspace'
419 0 : RETURN
420 : END IF
421 :
422 : ! the C pointer to the actual data
423 0 : buffer = C_LOC(attr_data)
424 :
425 : ! set the type of data
426 0 : type_id = h5t_native_integer
427 :
428 : ! create the attribute
429 0 : CALL h5acreate_f(loc_id, attr_name, type_id, space_id, attr_id, error)
430 0 : IF (error < 0) THEN
431 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
432 0 : ' ERROR: failed to create HDF5 attribute'
433 0 : RETURN
434 : END IF
435 :
436 : ! write the string attribute to file
437 0 : CALL h5awrite_f(attr_id, type_id, buffer, error)
438 0 : IF (error < 0) THEN
439 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
440 0 : ' ERROR: failed to write HDF5 attribute'
441 0 : RETURN
442 : END IF
443 :
444 : ! close attribute
445 0 : CALL h5aclose_f(attr_id, error)
446 0 : IF (error < 0) THEN
447 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
448 0 : ' ERROR: failed to close HDF5 attribute'
449 0 : RETURN
450 : END IF
451 :
452 : ! close dataspace
453 0 : CALL h5sclose_f(space_id, error)
454 0 : IF (error < 0) THEN
455 : WRITE (UNIT=output_unit, FMT="(/,T5,A,/)") &
456 0 : ' ERROR: failed to close HDF5 dataspace'
457 0 : RETURN
458 : END IF
459 :
460 : END SUBROUTINE h5awrite_integer_scalar
461 :
462 : ! **************************************************************************************************
463 : !> \brief Write a (scalar) double precision attribute
464 : !> \param loc_id either file id or group id
465 : !> \param attr_name the name of the attribute
466 : !> \param attr_data the attribute data, i.e. the double to write
467 : ! **************************************************************************************************
468 0 : SUBROUTINE h5awrite_double_scalar(loc_id, attr_name, attr_data)
469 : INTEGER(KIND=hid_t), INTENT(IN) :: loc_id
470 : CHARACTER(LEN=*), INTENT(IN) :: attr_name
471 : REAL(KIND=dp), INTENT(IN), TARGET :: attr_data
472 :
473 : INTEGER :: error
474 : INTEGER(KIND=hid_t) :: attr_id, space_id, type_id
475 : TYPE(c_ptr) :: buffer
476 :
477 : ! create a scalar dataspace
478 0 : CALL h5screate_f(h5s_scalar_f, space_id, error)
479 0 : IF (error < 0) CPABORT('ERROR: failed to create HDF5 dataspace')
480 :
481 : ! the C pointer to the actual data
482 0 : buffer = C_LOC(attr_data)
483 :
484 : ! set the type of data
485 0 : type_id = h5t_native_double
486 :
487 : ! create the attribute
488 0 : CALL h5acreate_f(loc_id, attr_name, type_id, space_id, attr_id, error)
489 0 : IF (error < 0) CPABORT('ERROR: failed to create HDF5 attribute')
490 :
491 : ! write the string attribute to file
492 0 : CALL h5awrite_f(attr_id, type_id, buffer, error)
493 0 : IF (error < 0) CPABORT('ERROR: failed to write HDF5 attribute')
494 :
495 : ! close attribute
496 0 : CALL h5aclose_f(attr_id, error)
497 0 : IF (error < 0) CPABORT('ERROR: failed to close HDF5 attribute')
498 :
499 : ! close dataspace
500 0 : CALL h5sclose_f(space_id, error)
501 0 : IF (error < 0) CPABORT('ERROR: failed to close HDF5 dataspace')
502 :
503 0 : END SUBROUTINE h5awrite_double_scalar
504 :
505 : ! **************************************************************************************************
506 : !> \brief Write an array of fixed-length string attribute
507 : !> \param loc_id either file id or group id
508 : !> \param attr_name the name of the attribute
509 : !> \param attr_data the attribute data, i.e. the array of strings
510 : ! **************************************************************************************************
511 0 : SUBROUTINE h5awrite_string_simple(loc_id, attr_name, attr_data)
512 : INTEGER(KIND=hid_t), INTENT(IN) :: loc_id
513 : CHARACTER(LEN=*), INTENT(IN) :: attr_name
514 : CHARACTER(LEN=*), DIMENSION(:), INTENT(IN), TARGET :: attr_data
515 :
516 : INTEGER :: error
517 : INTEGER(KIND=hid_t) :: attr_id, space_id, type_id
518 : INTEGER(KIND=hsize_t), DIMENSION(2) :: dims
519 : TYPE(c_ptr) :: buffer
520 :
521 0 : dims(1) = LEN(attr_data(1), kind=hsize_t) ! length of a string entry
522 0 : dims(2) = SIZE(attr_data, kind=hsize_t) ! length of array of strings
523 :
524 : ! create a fixed-length string datatype
525 0 : CALL h5tcopy_f(h5t_c_s1, type_id, error)
526 0 : CALL h5tset_cset_f(type_id, h5t_cset_utf8_f, error)
527 0 : CALL h5tset_size_f(type_id, INT(dims(1), size_t), error)
528 :
529 : ! create a simple dataspace
530 0 : CALL h5screate_simple_f(1, dims(2:2), space_id, error)
531 0 : IF (error < 0) CPABORT('ERROR: failed to create HDF5 dataspace')
532 :
533 : ! create the atrtibute
534 0 : CALL h5acreate_f(loc_id, attr_name, type_id, space_id, attr_id, error)
535 0 : IF (error < 0) CPABORT('ERROR: failed to create HDF5 attribute')
536 :
537 : ! the actual pointer to be passed
538 0 : buffer = C_LOC(attr_data(1))
539 :
540 : ! write the string array attribute to file
541 0 : CALL h5awrite_f(attr_id, type_id, buffer, error)
542 0 : IF (error < 0) CPABORT('ERROR: failed to write HDF5 attribute')
543 :
544 : ! close attribute
545 0 : CALL h5aclose_f(attr_id, error)
546 0 : IF (error < 0) CPABORT('ERROR: failed to close HDF5 attribute')
547 :
548 : ! close dataspace
549 0 : CALL h5sclose_f(space_id, error)
550 0 : IF (error < 0) CPABORT('ERROR: failed to close HDF5 dataspace')
551 :
552 : ! close datatype
553 0 : CALL h5tclose_f(type_id, error)
554 0 : IF (error < 0) CPABORT('ERROR: failed to close HDF5 datatype')
555 :
556 0 : END SUBROUTINE h5awrite_string_simple
557 :
558 : ! **************************************************************************************************
559 : !> \brief Write an array of doubles attribute
560 : !> \param loc_id either file id or group id
561 : !> \param attr_name the name of the attribute
562 : !> \param attr_data the attribute data, i.e. the array of doubles
563 : ! **************************************************************************************************
564 0 : SUBROUTINE h5awrite_double_simple(loc_id, attr_name, attr_data)
565 : INTEGER(KIND=hid_t), INTENT(IN) :: loc_id
566 : CHARACTER(LEN=*), INTENT(IN) :: attr_name
567 : REAL(KIND=dp), DIMENSION(:), INTENT(IN), TARGET :: attr_data
568 :
569 : INTEGER :: error
570 : INTEGER(KIND=hid_t) :: attr_id, space_id, type_id
571 : INTEGER(KIND=hsize_t), DIMENSION(1) :: dims
572 : TYPE(c_ptr) :: buffer
573 :
574 0 : dims(1) = SIZE(attr_data, kind=hsize_t) ! length of array of strings
575 :
576 : ! set the type of data
577 0 : type_id = h5t_native_double
578 :
579 : ! create a simple dataspace
580 0 : CALL h5screate_simple_f(1, dims, space_id, error)
581 0 : IF (error < 0) CPABORT('ERROR: failed to create HDF5 dataspace')
582 :
583 : ! create the atrtibute
584 0 : CALL h5acreate_f(loc_id, attr_name, type_id, space_id, attr_id, error)
585 0 : IF (error < 0) CPABORT('ERROR: failed to create HDF5 attribute')
586 :
587 : ! the actual pointer to be passed
588 0 : buffer = C_LOC(attr_data(1))
589 :
590 : ! write the string array attribute to file
591 0 : CALL h5awrite_f(attr_id, type_id, buffer, error)
592 0 : IF (error < 0) CPABORT('ERROR: failed to write HDF5 attribute')
593 :
594 : ! close attribute
595 0 : CALL h5aclose_f(attr_id, error)
596 0 : IF (error < 0) CPABORT('ERROR: failed to close HDF5 attribute')
597 :
598 : ! close dataspace
599 0 : CALL h5sclose_f(space_id, error)
600 0 : IF (error < 0) CPABORT('ERROR: failed to close HDF5 dataspace')
601 :
602 0 : END SUBROUTINE h5awrite_double_simple
603 :
604 : ! **************************************************************************************************
605 : !> \brief Write an array of integers attribute
606 : !> \param loc_id either file id or group id
607 : !> \param attr_name the name of the attribute
608 : !> \param attr_data the attribute data, i.e. the array of integers
609 : ! **************************************************************************************************
610 0 : SUBROUTINE h5awrite_integer_simple(loc_id, attr_name, attr_data)
611 : INTEGER(KIND=hid_t), INTENT(IN) :: loc_id
612 : CHARACTER(LEN=*), INTENT(IN) :: attr_name
613 : INTEGER, DIMENSION(:), INTENT(IN), TARGET :: attr_data
614 :
615 : INTEGER :: error
616 : INTEGER(KIND=hid_t) :: attr_id, space_id, type_id
617 : INTEGER(KIND=hsize_t), DIMENSION(1) :: dims
618 : TYPE(c_ptr) :: buffer
619 :
620 0 : dims(1) = SIZE(attr_data, kind=hsize_t) ! length of array of strings
621 :
622 : ! set the type of data
623 0 : type_id = h5t_native_integer
624 :
625 : ! create a simple dataspace
626 0 : CALL h5screate_simple_f(1, dims, space_id, error)
627 0 : IF (error < 0) CPABORT('ERROR: failed to create HDF5 dataspace')
628 :
629 : ! create the atrtibute
630 0 : CALL h5acreate_f(loc_id, attr_name, type_id, space_id, attr_id, error)
631 0 : IF (error < 0) CPABORT('ERROR: failed to create HDF5 attribute')
632 :
633 : ! the actual pointer to be passed
634 0 : buffer = C_LOC(attr_data(1))
635 :
636 : ! write the string array attribute to file
637 0 : CALL h5awrite_f(attr_id, type_id, buffer, error)
638 0 : IF (error < 0) CPABORT('ERROR: failed to write HDF5 attribute')
639 :
640 : ! close attribute
641 0 : CALL h5aclose_f(attr_id, error)
642 0 : IF (error < 0) CPABORT('ERROR: failed to close HDF5 attribute')
643 :
644 : ! close dataspace
645 0 : CALL h5sclose_f(space_id, error)
646 0 : IF (error < 0) CPABORT('ERROR: failed to close HDF5 dataspace')
647 :
648 0 : END SUBROUTINE h5awrite_integer_simple
649 :
650 : ! **************************************************************************************************
651 : !> \brief Write a dataset containing an array of doubles
652 : !> \param loc_id either file id or group id
653 : !> \param dset_name the name of the dataset
654 : !> \param dset_data the dataset data, i.e. the array of doubles
655 : ! **************************************************************************************************
656 0 : SUBROUTINE h5dwrite_double_simple(loc_id, dset_name, dset_data)
657 : INTEGER(KIND=hid_t), INTENT(IN) :: loc_id
658 : CHARACTER(LEN=*), INTENT(IN) :: dset_name
659 : REAL(KIND=dp), DIMENSION(:), INTENT(IN), TARGET :: dset_data
660 :
661 : INTEGER :: error
662 : INTEGER(KIND=hid_t) :: dset_id, space_id, type_id
663 : INTEGER(KIND=hsize_t), DIMENSION(1) :: dims
664 : TYPE(c_ptr) :: buffer
665 :
666 0 : dims(1) = SIZE(dset_data, kind=hsize_t) ! length of array
667 :
668 : ! set the type of data
669 0 : type_id = h5t_native_double
670 :
671 : ! create a simple dataspace
672 0 : CALL h5screate_simple_f(1, dims, space_id, error)
673 0 : IF (error < 0) CPABORT('ERROR: failed to create HDF5 dataspace')
674 :
675 : ! create the dataset
676 0 : CALL h5dcreate_f(loc_id, dset_name, type_id, space_id, dset_id, error)
677 0 : IF (error < 0) CPABORT('ERROR: failed to create HDF5 dataset')
678 :
679 : ! the actual pointer to be passed
680 0 : buffer = C_LOC(dset_data(1))
681 :
682 : ! write the string array attribute to file
683 0 : CALL h5dwrite_f(dset_id, type_id, buffer, error)
684 0 : IF (error < 0) CPABORT('ERROR: failed to write HDF5 dataset')
685 :
686 : ! close dataset
687 0 : CALL h5dclose_f(dset_id, error)
688 0 : IF (error < 0) CPABORT('ERROR: failed to close HDF5 dataset')
689 :
690 : ! close dataspace
691 0 : CALL h5sclose_f(space_id, error)
692 0 : IF (error < 0) CPABORT('ERROR: failed to close HDF5 dataspace')
693 :
694 0 : END SUBROUTINE h5dwrite_double_simple
695 :
696 : ! **************************************************************************************************
697 : !> \brief Read a dataset containing an array of doubles
698 : !> \param loc_id either file id or group id
699 : !> \param dset_name the name of the dataset
700 : !> \param dset_data where the read dataset data will be written
701 : ! **************************************************************************************************
702 0 : SUBROUTINE h5dread_double_simple(loc_id, dset_name, dset_data)
703 : INTEGER(KIND=hid_t), INTENT(IN) :: loc_id
704 : CHARACTER(LEN=*), INTENT(IN) :: dset_name
705 : REAL(KIND=dp), DIMENSION(:), INTENT(OUT) :: dset_data
706 :
707 : INTEGER :: error
708 : INTEGER(KIND=hid_t) :: dset_id, npoints, space_id, type_id
709 : INTEGER(KIND=hsize_t), DIMENSION(1) :: dims
710 :
711 0 : dims(1) = SIZE(dset_data, kind=hsize_t) ! length of array
712 :
713 : ! set the type of data
714 0 : type_id = h5t_native_double
715 :
716 : ! open the dataset
717 0 : CALL h5dopen_f(loc_id, dset_name, dset_id, error)
718 0 : IF (error < 0) CPABORT('ERROR: failed to open HDF5 dataset')
719 :
720 : ! get information on the dataspace
721 0 : CALL h5dget_space_f(dset_id, space_id, error)
722 0 : IF (error < 0) CPABORT('ERROR: failed to fetch HDF5 dataspace info')
723 :
724 : ! get dataspace dims
725 0 : CALL h5sget_simple_extent_npoints_f(space_id, npoints, error)
726 0 : IF (error < 0) CPABORT('ERROR: failed to fetch HDF5 dataspace dimension')
727 :
728 : ! read the data
729 0 : CALL h5dread_f(dset_id, type_id, dset_data, dims, error)
730 0 : IF (error < 0) CPABORT('ERROR: failed to read HDF5 dataset')
731 :
732 : ! close dataset
733 0 : CALL h5dclose_f(dset_id, error)
734 0 : IF (error < 0) CPABORT('ERROR: failed to close HDF5 dataset')
735 :
736 : ! close dataspace
737 0 : CALL h5sclose_f(space_id, error)
738 0 : IF (error < 0) CPABORT('ERROR: failed to close HDF5 dataspace')
739 :
740 0 : END SUBROUTINE h5dread_double_simple
741 :
742 : ! **************************************************************************************************
743 : !> \brief Read an attribute containing a scalar double
744 : !> \param loc_id either file id or group id
745 : !> \param attr_name ...
746 : !> \param attr_data ...
747 : ! **************************************************************************************************
748 0 : SUBROUTINE h5aread_double_scalar(loc_id, attr_name, attr_data)
749 : INTEGER(KIND=hid_t), INTENT(IN) :: loc_id
750 : CHARACTER(LEN=*), INTENT(IN) :: attr_name
751 : REAL(KIND=dp), INTENT(OUT), TARGET :: attr_data
752 :
753 : INTEGER :: error
754 : INTEGER(KIND=hid_t) :: attr_id, type_id
755 : TYPE(c_ptr) :: buffer
756 :
757 : ! set the type of data
758 0 : type_id = h5t_native_double
759 :
760 : ! open the attribute
761 0 : CALL h5aopen_f(loc_id, attr_name, attr_id, error)
762 0 : IF (error < 0) CPABORT('ERROR: failed to open HDF5 attribute')
763 :
764 0 : buffer = C_LOC(attr_data)
765 : ! read the data
766 0 : CALL h5aread_f(attr_id, type_id, buffer, error)
767 0 : IF (error < 0) CPABORT('ERROR: failed to read HDF5 attribute')
768 :
769 : ! close the attribute
770 0 : CALL h5aclose_f(attr_id, error)
771 0 : IF (error < 0) CPABORT('ERROR: failed to close HDF5 attribute')
772 :
773 0 : END SUBROUTINE h5aread_double_scalar
774 :
775 : #endif
776 :
777 : END MODULE hdf5_wrapper
|