29 #include "sll_assert.h"
30 #include "sll_working_precision.h"
53 h5oget_info_by_name_f, &
61 h5t_native_character, &
87 integer(hid_t),
private :: file_id
169 character(len=*),
intent(in) :: filename
171 integer,
intent(out) :: error
174 sll_assert_always(error == 0)
175 call h5fcreate_f(filename, h5f_acc_trunc_f, handle%file_id, error)
176 sll_assert_always(error == 0)
187 character(len=*),
intent(in) :: filename
189 integer,
intent(out) :: error
192 sll_assert_always(error == 0)
193 call h5fopen_f(filename, h5f_acc_rdwr_f, handle%file_id, error)
194 sll_assert_always(error == 0)
205 integer,
intent(out) :: error
207 call h5fclose_f(handle%file_id, error)
208 sll_assert_always(error == 0)
209 call h5close_f(error)
210 sll_assert_always(error == 0)
217 integer,
parameter :: rank = 1
219 integer(i32),
intent(in) :: array(:)
220 character(len=*),
intent(in) :: dsetname
221 integer,
intent(out) :: error
223 #define DATATYPE H5T_NATIVE_INTEGER
224 #include "sll_k_hdf5_ser_write_array.F90"
233 integer,
parameter :: rank = 2
235 integer(i32),
intent(in) :: array(:, :)
236 character(len=*),
intent(in) :: dsetname
237 integer,
intent(out) :: error
239 #define DATATYPE H5T_NATIVE_INTEGER
240 #include "sll_k_hdf5_ser_write_array.F90"
249 integer,
parameter :: rank = 3
251 integer(i32),
intent(in) :: array(:, :, :)
252 character(len=*),
intent(in) :: dsetname
253 integer,
intent(out) :: error
255 #define DATATYPE H5T_NATIVE_INTEGER
256 #include "sll_k_hdf5_ser_write_array.F90"
265 integer,
parameter :: rank = 1
267 real(f64),
intent(in) :: array(:)
268 character(len=*),
intent(in) :: dsetname
269 integer,
intent(out) :: error
271 #define DATATYPE H5T_NATIVE_DOUBLE
272 #include "sll_k_hdf5_ser_write_array.F90"
281 integer,
parameter :: rank = 2
283 real(f64),
intent(in) :: array(:, :)
284 character(len=*),
intent(in) :: dsetname
285 integer,
intent(out) :: error
287 #define DATATYPE H5T_NATIVE_DOUBLE
288 #include "sll_k_hdf5_ser_write_array.F90"
297 integer,
parameter :: rank = 3
299 real(f64),
intent(in) :: array(:, :, :)
300 character(len=*),
intent(in) :: dsetname
301 integer,
intent(out) :: error
303 #define DATATYPE H5T_NATIVE_DOUBLE
304 #include "sll_k_hdf5_ser_write_array.F90"
313 integer,
parameter :: rank = 1
315 real(f64),
intent(out) :: array(:)
316 character(len=*),
intent(in) :: dsetname
317 integer,
intent(out) :: error
319 #define DATATYPE H5T_NATIVE_DOUBLE
320 #include "sll_k_hdf5_ser_read_array.F90"
329 integer,
parameter :: rank = 2
331 real(f64),
intent(out) :: array(:, :)
332 character(len=*),
intent(in) :: dsetname
333 integer,
intent(out) :: error
335 #define DATATYPE H5T_NATIVE_DOUBLE
336 #include "sll_k_hdf5_ser_read_array.F90"
345 integer,
parameter :: rank = 3
347 real(f64),
intent(out) :: array(:, :, :)
348 character(len=*),
intent(in) :: dsetname
349 integer,
intent(out) :: error
351 #define DATATYPE H5T_NATIVE_DOUBLE
352 #include "sll_k_hdf5_ser_read_array.F90"
367 character(len=*),
intent(in) :: string
368 character(len=*),
intent(in) :: dsetname
369 integer,
intent(out) :: error
371 integer(hsize_t) :: array_dim(1)
372 integer(hid_t) :: dataset_id
373 integer(hid_t) :: dataspace_id
375 array_dim = int(len(string), hsize_t)
378 call h5screate_simple_f(1, array_dim, dataspace_id, error)
379 sll_assert_always(error == 0)
382 call h5dcreate_f(handle%file_id, dsetname, h5t_native_character, &
383 dataspace_id, dataset_id, error)
384 sll_assert_always(error == 0)
387 call h5dwrite_f(dataset_id, h5t_native_character, string, array_dim, error)
388 sll_assert_always(error == 0)
391 call h5sclose_f(dataspace_id, error)
392 sll_assert_always(error == 0)
395 call h5dclose_f(dataset_id, error)
396 sll_assert_always(error == 0)
410 character(len=*),
intent(in) :: filename
411 character(len=*),
intent(in) :: dsetname
412 integer,
intent(out) :: error
414 integer :: fileid, filesize
415 character(len=:),
allocatable :: content
417 open (newunit=fileid, file=filename, &
418 form=
'UNFORMATTED', access=
'STREAM', iostat=error)
421 inquire (file=filename, size=filesize)
422 if (filesize > 0)
then
425 allocate (
character(len=filesize) :: content)
428 read (fileid, pos=1, iostat=error) content
432 read (fileid, pos=filesize + 1, iostat=error) content
433 if (.not. is_iostat_end(error))
then
434 write (*, *)
"ERROR: file ", filename,
"was not completely read."
437 write (*, *)
'ERROR reading file: ', filename
441 write (*, *)
'Error getting size of file: ', filename
454 character(len=*),
intent(in) :: objpath
455 character(len=*),
intent(in) :: attrname
456 real(f64),
intent(in) :: attrvalue
457 integer,
intent(out) :: error
459 #define DATATYPE H5T_NATIVE_DOUBLE
460 #include "sll_k_hdf5_ser_write_attribute.F90"
470 character(len=*),
intent(in) :: objpath
471 character(len=*),
intent(in) :: attrname
472 integer,
intent(in) :: attrvalue
473 integer,
intent(out) :: error
475 #define DATATYPE H5T_NATIVE_INTEGER
476 #include "sll_k_hdf5_ser_write_attribute.F90"
486 character(len=*),
intent(in) :: objpath
487 character(len=*),
intent(in) :: attrname
488 real(f64),
intent(out) :: attrvalue
489 integer,
intent(out) :: error
491 #define DATATYPE H5T_NATIVE_DOUBLE
492 #include "sll_k_hdf5_ser_read_attribute.F90"
502 character(len=*),
intent(in) :: objpath
503 character(len=*),
intent(in) :: attrname
504 integer,
intent(out) :: attrvalue
505 integer,
intent(out) :: error
507 #define DATATYPE H5T_NATIVE_INTEGER
508 #include "sll_k_hdf5_ser_read_attribute.F90"
Read nD array of double precision floats or integers from HDF5 file.
Read pre-existing named scalar attribute (double precision float or integer) from HDF5 object (group ...
Write nD array of double precision floats or integers into HDF5 file.
Attach new named scalar attribute (double precision float or integer) to HDF5 object (group or datase...
Implements the functions to write hdf5 file to store heavy data.
subroutine sll_hdf5_ser_write_dble_array_2d(handle, array, dsetname, error)
Write 2D array of float64 into HDF5 file.
subroutine s_hdf5_ser_write_attribute_dble(handle, objpath, attrname, attrvalue, error)
Attach named float64 attribute to existing group or dataset.
subroutine sll_hdf5_ser_write_char_array(handle, string, dsetname, error)
Write Fortran string to HDF5 file as 1D array of characters.
integer, parameter h5o_type_group
subroutine sll_hdf5_ser_read_dble_array_1d(handle, array, dsetname, error)
Read 1D array of float64 from HDF5 file.
subroutine s_hdf5_ser_read_attribute_int(handle, objpath, attrname, attrvalue, error)
Read named integer attribute from existing group or dataset.
subroutine s_hdf5_ser_read_attribute_dble(handle, objpath, attrname, attrvalue, error)
Read named float64 attribute from existing group or dataset.
subroutine sll_hdf5_ser_write_int_array_1d(handle, array, dsetname, error)
Write 1D array of int32 into HDF5 file.
subroutine, public sll_s_hdf5_ser_file_create(filename, handle, error)
Create new HDF5 file.
subroutine sll_hdf5_ser_write_dble_array_1d(handle, array, dsetname, error)
Write 1D array of float64 into HDF5 file.
subroutine sll_hdf5_ser_write_int_array_3d(handle, array, dsetname, error)
Write 3D array of int32 into HDF5 file.
subroutine sll_hdf5_ser_write_int_array_2d(handle, array, dsetname, error)
Write 2D array of int32 into HDF5 file.
subroutine sll_hdf5_ser_read_dble_array_2d(handle, array, dsetname, error)
Read 2D array of float64 from HDF5 file.
subroutine sll_hdf5_ser_read_dble_array_3d(handle, array, dsetname, error)
Read 3D array of float64 from HDF5 file.
subroutine s_hdf5_ser_write_attribute_int(handle, objpath, attrname, attrvalue, error)
Attach named integer attribute to existing group or dataset.
integer, parameter h5o_type_dataset
subroutine, public sll_s_hdf5_ser_file_close(handle, error)
Close existing HDF5 file.
subroutine, public sll_s_hdf5_ser_file_open(filename, handle, error)
Open existing HDF5 file.
subroutine, public sll_s_hdf5_ser_write_file(handle, filename, dsetname, error)
Read complete file to string and store it into HDF5 data structure.
subroutine sll_hdf5_ser_write_dble_array_3d(handle, array, dsetname, error)
Write 3D array of float64 into HDF5 file.
Opaque object around HDF5 file id.