Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_hdf5_io_serial.F90
Go to the documentation of this file.
1 !**************************************************************
2 ! Copyright INRIA
3 !
4 ! This code SeLaLib (for Semi-Lagrangian-Library)
5 ! is a parallel library for simulating the plasma turbulence
6 ! in a tokamak.
7 !
8 ! This software is governed by the CeCILL-B license
9 ! under French law and abiding by the rules of distribution
10 ! of free software. You can use, modify and redistribute
11 ! the software under the terms of the CeCILL-B license as
12 ! circulated by CEA, CNRS and INRIA at the following URL
13 ! "http://www.cecill.info".
14 !**************************************************************
15 
26 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
27 #ifndef NOHDF5
28 
29 #include "sll_assert.h"
30 #include "sll_working_precision.h"
31 
32  use hdf5, only: &
33  h5acreate_f, &
34  h5aclose_f, &
35  h5aopen_f, &
36  h5aread_f, &
37  h5awrite_f, &
38  h5close_f, &
39  h5dclose_f, &
40  h5dcreate_f, &
41  h5dopen_f, &
42  h5dread_f, &
43  h5dwrite_f, &
44  h5f_acc_rdonly_f, &
45  h5f_acc_rdwr_f, &
46  h5f_acc_trunc_f, &
47  h5fclose_f, &
48  h5fcreate_f, &
49  h5fopen_f, &
50  h5gclose_f, &
51  h5gopen_f, &
52  h5o_info_t, &
53  h5oget_info_by_name_f, &
54  h5open_f, &
55  h5s_scalar_f, &
56  h5sclose_f, &
57  h5screate_f, &
58  h5screate_simple_f, &
59  h5t_native_double, &
60  h5t_native_integer, &
61  h5t_native_character, &
62  hid_t, &
63  hsize_t
64 
65  implicit none
66 
67  public :: &
77 
78  private
79 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
80 
81  ! HDF5 workaround: C enumeration does not seem to be available in Fortran
82  integer, parameter :: h5o_type_group = 0
83  integer, parameter :: h5o_type_dataset = 1
84 
87  integer(hid_t), private :: file_id
88  end type
89 
90  !-----------------------------------------------------------------------------
98  !-----------------------------------------------------------------------------
100  module procedure sll_hdf5_ser_write_dble_array_1d
101  module procedure sll_hdf5_ser_write_dble_array_2d
102  module procedure sll_hdf5_ser_write_dble_array_3d
103  module procedure sll_hdf5_ser_write_int_array_1d
104  module procedure sll_hdf5_ser_write_int_array_2d
105  module procedure sll_hdf5_ser_write_int_array_3d
106  end interface
107 
108  !-----------------------------------------------------------------------------
116  !-----------------------------------------------------------------------------
118  module procedure sll_hdf5_ser_read_dble_array_1d
119  module procedure sll_hdf5_ser_read_dble_array_2d
120  module procedure sll_hdf5_ser_read_dble_array_3d
121  end interface
122 
123  !-----------------------------------------------------------------------------
133  !-----------------------------------------------------------------------------
135  module procedure s_hdf5_ser_write_attribute_dble
136  module procedure s_hdf5_ser_write_attribute_int
137  end interface
138 
139  !-----------------------------------------------------------------------------
149  !-----------------------------------------------------------------------------
151  module procedure s_hdf5_ser_read_attribute_dble
152  module procedure s_hdf5_ser_read_attribute_int
153  end interface
154 
155 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
156 contains
157 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
158 
159  !-----------------------------------------------------------------------------
167  !-----------------------------------------------------------------------------
168  subroutine sll_s_hdf5_ser_file_create(filename, handle, error)
169  character(len=*), intent(in) :: filename
170  type(sll_t_hdf5_ser_handle), intent(out) :: handle
171  integer, intent(out) :: error
172 
173  call h5open_f(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)
177  end subroutine sll_s_hdf5_ser_file_create
178 
179  !-----------------------------------------------------------------------------
185  !-----------------------------------------------------------------------------
186  subroutine sll_s_hdf5_ser_file_open(filename, handle, error)
187  character(len=*), intent(in) :: filename
188  type(sll_t_hdf5_ser_handle), intent(out) :: handle
189  integer, intent(out) :: error
190 
191  call h5open_f(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)
195  end subroutine sll_s_hdf5_ser_file_open
196 
197  !-----------------------------------------------------------------------------
202  !-----------------------------------------------------------------------------
203  subroutine sll_s_hdf5_ser_file_close(handle, error)
204  type(sll_t_hdf5_ser_handle), intent(in) :: handle
205  integer, intent(out) :: error
206 
207  call h5fclose_f(handle%file_id, error) ! Close property list and file
208  sll_assert_always(error == 0)
209  call h5close_f(error) ! Close FORTRAN interface
210  sll_assert_always(error == 0)
211  end subroutine sll_s_hdf5_ser_file_close
212 
213  !-----------------------------------------------------------------------------
215  !-----------------------------------------------------------------------------
216  subroutine sll_hdf5_ser_write_int_array_1d(handle, array, dsetname, error)
217  integer, parameter :: rank = 1
218  type(sll_t_hdf5_ser_handle), intent(in) :: handle
219  integer(i32), intent(in) :: array(:)
220  character(len=*), intent(in) :: dsetname
221  integer, intent(out) :: error
222 
223 #define DATATYPE H5T_NATIVE_INTEGER
224 #include "sll_k_hdf5_ser_write_array.F90"
225 #undef DATATYPE
226 
227  end subroutine sll_hdf5_ser_write_int_array_1d
228 
229  !-----------------------------------------------------------------------------
231  !-----------------------------------------------------------------------------
232  subroutine sll_hdf5_ser_write_int_array_2d(handle, array, dsetname, error)
233  integer, parameter :: rank = 2
234  type(sll_t_hdf5_ser_handle), intent(in) :: handle
235  integer(i32), intent(in) :: array(:, :)
236  character(len=*), intent(in) :: dsetname
237  integer, intent(out) :: error
238 
239 #define DATATYPE H5T_NATIVE_INTEGER
240 #include "sll_k_hdf5_ser_write_array.F90"
241 #undef DATATYPE
242 
243  end subroutine sll_hdf5_ser_write_int_array_2d
244 
245  !-----------------------------------------------------------------------------
247  !-----------------------------------------------------------------------------
248  subroutine sll_hdf5_ser_write_int_array_3d(handle, array, dsetname, error)
249  integer, parameter :: rank = 3
250  type(sll_t_hdf5_ser_handle), intent(in) :: handle
251  integer(i32), intent(in) :: array(:, :, :)
252  character(len=*), intent(in) :: dsetname
253  integer, intent(out) :: error
254 
255 #define DATATYPE H5T_NATIVE_INTEGER
256 #include "sll_k_hdf5_ser_write_array.F90"
257 #undef DATATYPE
258 
259  end subroutine sll_hdf5_ser_write_int_array_3d
260 
261  !-----------------------------------------------------------------------------
263  !-----------------------------------------------------------------------------
264  subroutine sll_hdf5_ser_write_dble_array_1d(handle, array, dsetname, error)
265  integer, parameter :: rank = 1
266  type(sll_t_hdf5_ser_handle), intent(in) :: handle
267  real(f64), intent(in) :: array(:)
268  character(len=*), intent(in) :: dsetname
269  integer, intent(out) :: error
270 
271 #define DATATYPE H5T_NATIVE_DOUBLE
272 #include "sll_k_hdf5_ser_write_array.F90"
273 #undef DATATYPE
274 
275  end subroutine sll_hdf5_ser_write_dble_array_1d
276 
277  !-----------------------------------------------------------------------------
279  !-----------------------------------------------------------------------------
280  subroutine sll_hdf5_ser_write_dble_array_2d(handle, array, dsetname, error)
281  integer, parameter :: rank = 2
282  type(sll_t_hdf5_ser_handle), intent(in) :: handle
283  real(f64), intent(in) :: array(:, :)
284  character(len=*), intent(in) :: dsetname
285  integer, intent(out) :: error
286 
287 #define DATATYPE H5T_NATIVE_DOUBLE
288 #include "sll_k_hdf5_ser_write_array.F90"
289 #undef DATATYPE
290 
291  end subroutine sll_hdf5_ser_write_dble_array_2d
292 
293  !-----------------------------------------------------------------------------
295  !-----------------------------------------------------------------------------
296  subroutine sll_hdf5_ser_write_dble_array_3d(handle, array, dsetname, error)
297  integer, parameter :: rank = 3
298  type(sll_t_hdf5_ser_handle), intent(in) :: handle
299  real(f64), intent(in) :: array(:, :, :)
300  character(len=*), intent(in) :: dsetname
301  integer, intent(out) :: error
302 
303 #define DATATYPE H5T_NATIVE_DOUBLE
304 #include "sll_k_hdf5_ser_write_array.F90"
305 #undef DATATYPE
306 
307  end subroutine sll_hdf5_ser_write_dble_array_3d
308 
309  !-----------------------------------------------------------------------------
311  !-----------------------------------------------------------------------------
312  subroutine sll_hdf5_ser_read_dble_array_1d(handle, array, dsetname, error)
313  integer, parameter :: rank = 1
314  type(sll_t_hdf5_ser_handle), intent(in) :: handle
315  real(f64), intent(out) :: array(:)
316  character(len=*), intent(in) :: dsetname
317  integer, intent(out) :: error
318 
319 #define DATATYPE H5T_NATIVE_DOUBLE
320 #include "sll_k_hdf5_ser_read_array.F90"
321 #undef DATATYPE
322 
323  end subroutine sll_hdf5_ser_read_dble_array_1d
324 
325  !-----------------------------------------------------------------------------
327  !-----------------------------------------------------------------------------
328  subroutine sll_hdf5_ser_read_dble_array_2d(handle, array, dsetname, error)
329  integer, parameter :: rank = 2
330  type(sll_t_hdf5_ser_handle), intent(in) :: handle
331  real(f64), intent(out) :: array(:, :)
332  character(len=*), intent(in) :: dsetname
333  integer, intent(out) :: error
334 
335 #define DATATYPE H5T_NATIVE_DOUBLE
336 #include "sll_k_hdf5_ser_read_array.F90"
337 #undef DATATYPE
338 
339  end subroutine sll_hdf5_ser_read_dble_array_2d
340 
341  !-----------------------------------------------------------------------------
343  !-----------------------------------------------------------------------------
344  subroutine sll_hdf5_ser_read_dble_array_3d(handle, array, dsetname, error)
345  integer, parameter :: rank = 3
346  type(sll_t_hdf5_ser_handle), intent(in) :: handle
347  real(f64), intent(out) :: array(:, :, :)
348  character(len=*), intent(in) :: dsetname
349  integer, intent(out) :: error
350 
351 #define DATATYPE H5T_NATIVE_DOUBLE
352 #include "sll_k_hdf5_ser_read_array.F90"
353 #undef DATATYPE
354 
355  end subroutine sll_hdf5_ser_read_dble_array_3d
356 
357  !-----------------------------------------------------------------------------
364  !-----------------------------------------------------------------------------
365  subroutine sll_hdf5_ser_write_char_array(handle, string, dsetname, error)
366  type(sll_t_hdf5_ser_handle), intent(in) :: handle
367  character(len=*), intent(in) :: string
368  character(len=*), intent(in) :: dsetname
369  integer, intent(out) :: error
370 
371  integer(hsize_t) :: array_dim(1)
372  integer(hid_t) :: dataset_id
373  integer(hid_t) :: dataspace_id
374 
375  array_dim = int(len(string), hsize_t)
376 
377  ! Create dataspace
378  call h5screate_simple_f(1, array_dim, dataspace_id, error)
379  sll_assert_always(error == 0)
380 
381  ! Create dataset
382  call h5dcreate_f(handle%file_id, dsetname, h5t_native_character, &
383  dataspace_id, dataset_id, error)
384  sll_assert_always(error == 0)
385 
386  ! Write dataset
387  call h5dwrite_f(dataset_id, h5t_native_character, string, array_dim, error)
388  sll_assert_always(error == 0)
389 
390  ! Close dataspace
391  call h5sclose_f(dataspace_id, error)
392  sll_assert_always(error == 0)
393 
394  ! Close dataset
395  call h5dclose_f(dataset_id, error)
396  sll_assert_always(error == 0)
397 
398  end subroutine sll_hdf5_ser_write_char_array
399 
400  !-----------------------------------------------------------------------------
407  !-----------------------------------------------------------------------------
408  subroutine sll_s_hdf5_ser_write_file(handle, filename, dsetname, error)
409  type(sll_t_hdf5_ser_handle), intent(in) :: handle
410  character(len=*), intent(in) :: filename
411  character(len=*), intent(in) :: dsetname
412  integer, intent(out) :: error
413 
414  integer :: fileid, filesize
415  character(len=:), allocatable :: content
416 
417  open (newunit=fileid, file=filename, &
418  form='UNFORMATTED', access='STREAM', iostat=error)
419 
420  !get size of file
421  inquire (file=filename, size=filesize)
422  if (filesize > 0) then
423 
424  ! Allocate memory
425  allocate (character(len=filesize) :: content)
426 
427  ! Read the whole file into one string
428  read (fileid, pos=1, iostat=error) content
429 
430  if (error == 0) then
431  !try to read more in case not everything was read
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."
435  end if
436  else
437  write (*, *) 'ERROR reading file: ', filename
438  end if
439 
440  else
441  write (*, *) 'Error getting size of file: ', filename
442  end if
443 
444  close (fileid)
445  call sll_hdf5_ser_write_char_array(handle, content, dsetname, error)
446 
447  end subroutine sll_s_hdf5_ser_write_file
448 
449  !-----------------------------------------------------------------------------
451  !-----------------------------------------------------------------------------
452  subroutine s_hdf5_ser_write_attribute_dble(handle, objpath, attrname, attrvalue, error)
453  type(sll_t_hdf5_ser_handle), intent(in) :: handle
454  character(len=*), intent(in) :: objpath
455  character(len=*), intent(in) :: attrname
456  real(f64), intent(in) :: attrvalue
457  integer, intent(out) :: error
458 
459 #define DATATYPE H5T_NATIVE_DOUBLE
460 #include "sll_k_hdf5_ser_write_attribute.F90"
461 #undef DATATYPE
462 
463  end subroutine s_hdf5_ser_write_attribute_dble
464 
465  !-----------------------------------------------------------------------------
467  !-----------------------------------------------------------------------------
468  subroutine s_hdf5_ser_write_attribute_int(handle, objpath, attrname, attrvalue, error)
469  type(sll_t_hdf5_ser_handle), intent(in) :: handle
470  character(len=*), intent(in) :: objpath
471  character(len=*), intent(in) :: attrname
472  integer, intent(in) :: attrvalue
473  integer, intent(out) :: error
474 
475 #define DATATYPE H5T_NATIVE_INTEGER
476 #include "sll_k_hdf5_ser_write_attribute.F90"
477 #undef DATATYPE
478 
479  end subroutine s_hdf5_ser_write_attribute_int
480 
481  !-----------------------------------------------------------------------------
483  !-----------------------------------------------------------------------------
484  subroutine s_hdf5_ser_read_attribute_dble(handle, objpath, attrname, attrvalue, error)
485  type(sll_t_hdf5_ser_handle), intent(in) :: handle
486  character(len=*), intent(in) :: objpath
487  character(len=*), intent(in) :: attrname
488  real(f64), intent(out) :: attrvalue
489  integer, intent(out) :: error
490 
491 #define DATATYPE H5T_NATIVE_DOUBLE
492 #include "sll_k_hdf5_ser_read_attribute.F90"
493 #undef DATATYPE
494 
495  end subroutine s_hdf5_ser_read_attribute_dble
496 
497  !-----------------------------------------------------------------------------
499  !-----------------------------------------------------------------------------
500  subroutine s_hdf5_ser_read_attribute_int(handle, objpath, attrname, attrvalue, error)
501  type(sll_t_hdf5_ser_handle), intent(in) :: handle
502  character(len=*), intent(in) :: objpath
503  character(len=*), intent(in) :: attrname
504  integer, intent(out) :: attrvalue
505  integer, intent(out) :: error
506 
507 #define DATATYPE H5T_NATIVE_INTEGER
508 #include "sll_k_hdf5_ser_read_attribute.F90"
509 #undef DATATYPE
510 
511  end subroutine s_hdf5_ser_read_attribute_int
512 
513  !-----------------------------------------------------------------------------
514 #endif
515 
516 end module sll_m_hdf5_io_serial
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.
    Report Typos and Errors