Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_hdf5_io_parallel.F90
Go to the documentation of this file.
1 !**************************************************************
2 ! Copyright INRIA
3 ! Authors :
4 ! Pierre Navaro
5 !
6 ! This code SeLaLib (for Semi-Lagrangian-Library)
7 ! is a parallel library for simulating the plasma turbulence
8 ! in a tokamak.
9 !
10 ! This software is governed by the CeCILL-B license
11 ! under French law and abiding by the rules of distribution
12 ! of free software. You can use, modify and redistribute
13 ! the software under the terms of the CeCILL-B license as
14 ! circulated by CEA, CNRS and INRIA at the following URL
15 ! "http://www.cecill.info".
16 !**************************************************************
17 
26 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
27 #ifndef NOHDF5
28 
29 #include "sll_assert.h"
30 #include "sll_working_precision.h"
31 
32  use hdf5, only: &
33  h5close_f, &
34  h5dclose_f, &
35  h5dcreate_f, &
36  h5dget_space_f, &
37  h5dopen_f, &
38  h5dread_f, &
39  h5dwrite_f, &
40  h5f_acc_rdonly_f, &
41  h5f_acc_trunc_f, &
42  h5fclose_f, &
43  h5fcreate_f, &
44  h5fd_mpio_collective_f, &
45  h5fopen_f, &
46  h5open_f, &
47  h5p_dataset_create_f, &
48  h5p_dataset_xfer_f, &
49  h5p_file_access_f, &
50  h5pclose_f, &
51  h5pcreate_f, &
52  h5pset_chunk_f, &
53  h5pset_dxpl_mpio_f, &
54  h5pset_fapl_mpio_f, &
55  h5s_select_and_f, &
56  h5s_select_set_f, &
57  h5sclose_f, &
58  h5screate_simple_f, &
59  h5sselect_hyperslab_f, &
60  h5t_native_double, &
61  hid_t, &
62  hsize_t, &
63  hssize_t
64 
65  use mpi, only: &
66  mpi_info_null
67 
68  implicit none
69 
70  public :: &
77 
78  private
79 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
80 
83  integer(hid_t), private :: file_id
84  end type
85 
86  !-----------------------------------------------------------------------------
101  !-----------------------------------------------------------------------------
103  module procedure sll_hdf5_par_write_dble_array_1d
104  module procedure sll_hdf5_par_write_dble_array_2d
105  module procedure sll_hdf5_par_write_dble_array_3d
106  module procedure sll_hdf5_par_write_dble_array_4d
107  module procedure sll_hdf5_par_write_dble_array_5d
108  module procedure sll_hdf5_par_write_dble_array_6d
109  end interface
110 
111  !-----------------------------------------------------------------------------
125  !-----------------------------------------------------------------------------
127  module procedure sll_hdf5_par_read_dble_array_1d
128  module procedure sll_hdf5_par_read_dble_array_2d
129  module procedure sll_hdf5_par_read_dble_array_3d
130  module procedure sll_hdf5_par_read_dble_array_4d
131  module procedure sll_hdf5_par_read_dble_array_5d
132  module procedure sll_hdf5_par_read_dble_array_6d
133  end interface sll_o_hdf5_par_read_array
134 
135 contains
136 
137  !-----------------------------------------------------------------------------
146  !-----------------------------------------------------------------------------
147  subroutine sll_s_hdf5_par_file_create(filename, comm, handle, error)
148  character(len=*), intent(in) :: filename
149  integer, intent(in) :: comm
150  type(sll_t_hdf5_par_handle), intent(out) :: handle
151  integer, intent(out) :: error
152 
153  integer(hid_t) :: plist_id
154  integer :: info
155  info = mpi_info_null
156 
157  call h5open_f(error)
158  sll_assert_always(error == 0)
159  call h5pcreate_f(h5p_file_access_f, plist_id, error)
160  sll_assert_always(error == 0)
161  call h5pset_fapl_mpio_f(plist_id, comm, info, error)
162  sll_assert_always(error == 0)
163  call h5fcreate_f(filename, h5f_acc_trunc_f, handle%file_id, error, access_prp=plist_id)
164  sll_assert_always(error == 0)
165  call h5pclose_f(plist_id, error)
166  sll_assert_always(error == 0)
167 
168  end subroutine sll_s_hdf5_par_file_create
169 
170  !-----------------------------------------------------------------------------
179  !-----------------------------------------------------------------------------
180  subroutine sll_s_hdf5_par_file_open(filename, comm, handle, error)
181  character(len=*), intent(in) :: filename
182  integer, intent(in) :: comm
183  type(sll_t_hdf5_par_handle), intent(out) :: handle
184  integer, intent(out) :: error
185 
186  integer(hid_t) :: plist_id
187  integer :: info
188  info = mpi_info_null
189 
190  call h5open_f(error)
191  sll_assert_always(error == 0)
192  call h5pcreate_f(h5p_file_access_f, plist_id, error)
193  sll_assert_always(error == 0)
194  call h5pset_fapl_mpio_f(plist_id, comm, info, error)
195  sll_assert_always(error == 0)
196  call h5fopen_f(filename, h5f_acc_rdonly_f, handle%file_id, error, access_prp=plist_id)
197  sll_assert_always(error == 0)
198  call h5pclose_f(plist_id, error)
199  sll_assert_always(error == 0)
200 
201  end subroutine sll_s_hdf5_par_file_open
202 
203  !-----------------------------------------------------------------------------
208  !-----------------------------------------------------------------------------
209  subroutine sll_s_hdf5_par_file_close(handle, error)
210  type(sll_t_hdf5_par_handle), intent(in) :: handle
211  integer, intent(out) :: error
212 
213  call h5fclose_f(handle%file_id, error) ! Close property list and file
214  sll_assert_always(error == 0)
215  call h5close_f(error) ! Close FORTRAN interface
216  sll_assert_always(error == 0)
217  end subroutine sll_s_hdf5_par_file_close
218 
219  !-----------------------------------------------------------------------------
221  !-----------------------------------------------------------------------------
223  handle, global_size, offset, array, dsetname, error, chunk_dims)
224  integer, parameter :: dspace_dims = 1
225  type(sll_t_hdf5_par_handle), intent(in) :: handle
226  integer(i64), intent(in) :: global_size(dspace_dims)
227  integer(i64), intent(in) :: offset(dspace_dims)
228  real(f64), intent(in) :: array(:)
229  character(len=*), intent(in) :: dsetname
230  integer, intent(out) :: error
231  integer(i64), optional, intent(in) :: chunk_dims(dspace_dims)
232 
233 #define DATATYPE H5T_NATIVE_DOUBLE
234 #include "sll_k_hdf5_par_write_array.F90"
235 #undef DATATYPE
236 
237  end subroutine sll_hdf5_par_write_dble_array_1d
238 
239  !-----------------------------------------------------------------------------
241  !-----------------------------------------------------------------------------
243  handle, global_size, offset, array, dsetname, error, chunk_dims)
244  integer, parameter :: dspace_dims = 2
245  type(sll_t_hdf5_par_handle), intent(in) :: handle
246  integer(i64), intent(in) :: global_size(dspace_dims)
247  integer(i64), intent(in) :: offset(dspace_dims)
248  real(f64), intent(in) :: array(:, :)
249  character(len=*), intent(in) :: dsetname
250  integer, intent(out) :: error
251  integer(i64), optional, intent(in) :: chunk_dims(dspace_dims)
252 
253 #define DATATYPE H5T_NATIVE_DOUBLE
254 #include "sll_k_hdf5_par_write_array.F90"
255 #undef DATATYPE
256 
257  end subroutine sll_hdf5_par_write_dble_array_2d
258 
259  !-----------------------------------------------------------------------------
261  !-----------------------------------------------------------------------------
263  handle, global_size, offset, array, dsetname, error, chunk_dims)
264  integer, parameter :: dspace_dims = 3
265  type(sll_t_hdf5_par_handle), intent(in) :: handle
266  integer(i64), intent(in) :: global_size(dspace_dims)
267  integer(i64), intent(in) :: offset(dspace_dims)
268  real(f64), intent(in) :: array(:, :, :)
269  character(len=*), intent(in) :: dsetname
270  integer, intent(out) :: error
271  integer(i64), optional, intent(in) :: chunk_dims(dspace_dims)
272 
273 #define DATATYPE H5T_NATIVE_DOUBLE
274 #include "sll_k_hdf5_par_write_array.F90"
275 #undef DATATYPE
276 
277  end subroutine sll_hdf5_par_write_dble_array_3d
278 
279  !-----------------------------------------------------------------------------
281  !-----------------------------------------------------------------------------
283  handle, global_size, offset, array, dsetname, error, chunk_dims)
284  integer, parameter :: dspace_dims = 4
285  type(sll_t_hdf5_par_handle), intent(in) :: handle
286  integer(i64), intent(in) :: global_size(dspace_dims)
287  integer(i64), intent(in) :: offset(dspace_dims)
288  real(f64), intent(in) :: array(:, :, :, :)
289  character(len=*), intent(in) :: dsetname
290  integer, intent(out) :: error
291  integer(i64), optional, intent(in) :: chunk_dims(dspace_dims)
292 
293 #define DATATYPE H5T_NATIVE_DOUBLE
294 #include "sll_k_hdf5_par_write_array.F90"
295 #undef DATATYPE
296 
297  end subroutine sll_hdf5_par_write_dble_array_4d
298 
299  !-----------------------------------------------------------------------------
301  !-----------------------------------------------------------------------------
303  handle, global_size, offset, array, dsetname, error, chunk_dims)
304  integer, parameter :: dspace_dims = 5
305  type(sll_t_hdf5_par_handle), intent(in) :: handle
306  integer(i64), intent(in) :: global_size(dspace_dims)
307  integer(i64), intent(in) :: offset(dspace_dims)
308  real(f64), intent(in) :: array(:, :, :, :, :)
309  character(len=*), intent(in) :: dsetname
310  integer, intent(out) :: error
311  integer(i64), optional, intent(in) :: chunk_dims(dspace_dims)
312 
313 #define DATATYPE H5T_NATIVE_DOUBLE
314 #include "sll_k_hdf5_par_write_array.F90"
315 #undef DATATYPE
316 
317  end subroutine sll_hdf5_par_write_dble_array_5d
318 
319  !-----------------------------------------------------------------------------
321  !-----------------------------------------------------------------------------
323  handle, global_size, offset, array, dsetname, error, chunk_dims)
324  integer, parameter :: dspace_dims = 6
325  type(sll_t_hdf5_par_handle), intent(in) :: handle
326  integer(i64), intent(in) :: global_size(dspace_dims)
327  integer(i64), intent(in) :: offset(dspace_dims)
328  real(f64), intent(in) :: array(:, :, :, :, :, :)
329  character(len=*), intent(in) :: dsetname
330  integer, intent(out) :: error
331  integer(i64), optional, intent(in) :: chunk_dims(dspace_dims)
332 
333 #define DATATYPE H5T_NATIVE_DOUBLE
334 #include "sll_k_hdf5_par_write_array.F90"
335 #undef DATATYPE
336 
337  end subroutine sll_hdf5_par_write_dble_array_6d
338 
339  !-----------------------------------------------------------------------------
341  !-----------------------------------------------------------------------------
343  handle, global_size, offset, array, dsetname, error)
344  integer, parameter :: dspace_dims = 1
345  type(sll_t_hdf5_par_handle), intent(in) :: handle
346  integer(i64), intent(in) :: global_size(dspace_dims)
347  integer(i64), intent(in) :: offset(dspace_dims)
348  real(f64), intent(out) :: array(:)
349  character(len=*), intent(in) :: dsetname
350  integer, intent(out) :: error
351 
352 #define DATATYPE H5T_NATIVE_DOUBLE
353 #include "sll_k_hdf5_par_read_array.F90"
354 #undef DATATYPE
355 
356  end subroutine sll_hdf5_par_read_dble_array_1d
357 
358  !-----------------------------------------------------------------------------
360  !-----------------------------------------------------------------------------
362  handle, global_size, offset, array, dsetname, error)
363  integer, parameter :: dspace_dims = 2
364  type(sll_t_hdf5_par_handle), intent(in) :: handle
365  integer(i64), intent(in) :: global_size(dspace_dims)
366  integer(i64), intent(in) :: offset(dspace_dims)
367  real(f64), intent(out) :: array(:, :)
368  character(len=*), intent(in) :: dsetname
369  integer, intent(out) :: error
370 
371 #define DATATYPE H5T_NATIVE_DOUBLE
372 #include "sll_k_hdf5_par_read_array.F90"
373 #undef DATATYPE
374 
375  end subroutine sll_hdf5_par_read_dble_array_2d
376 
377  !-----------------------------------------------------------------------------
379  !-----------------------------------------------------------------------------
381  handle, global_size, offset, array, dsetname, error)
382  integer, parameter :: dspace_dims = 3
383  type(sll_t_hdf5_par_handle), intent(in) :: handle
384  integer(i64), intent(in) :: global_size(dspace_dims)
385  integer(i64), intent(in) :: offset(dspace_dims)
386  real(f64), intent(out) :: array(:, :, :)
387  character(len=*), intent(in) :: dsetname
388  integer, intent(out) :: error
389 
390 #define DATATYPE H5T_NATIVE_DOUBLE
391 #include "sll_k_hdf5_par_read_array.F90"
392 #undef DATATYPE
393 
394  end subroutine sll_hdf5_par_read_dble_array_3d
395 
396  !-----------------------------------------------------------------------------
398  !-----------------------------------------------------------------------------
400  handle, global_size, offset, array, dsetname, error)
401  integer, parameter :: dspace_dims = 4
402  type(sll_t_hdf5_par_handle), intent(in) :: handle
403  integer(i64), intent(in) :: global_size(dspace_dims)
404  integer(i64), intent(in) :: offset(dspace_dims)
405  real(f64), intent(out) :: array(:, :, :, :)
406  character(len=*), intent(in) :: dsetname
407  integer, intent(out) :: error
408 
409 #define DATATYPE H5T_NATIVE_DOUBLE
410 #include "sll_k_hdf5_par_read_array.F90"
411 #undef DATATYPE
412 
413  end subroutine sll_hdf5_par_read_dble_array_4d
414 
415  !-----------------------------------------------------------------------------
417  !-----------------------------------------------------------------------------
419  handle, global_size, offset, array, dsetname, error)
420  integer, parameter :: dspace_dims = 5
421  type(sll_t_hdf5_par_handle), intent(in) :: handle
422  integer(i64), intent(in) :: global_size(dspace_dims)
423  integer(i64), intent(in) :: offset(dspace_dims)
424  real(f64), intent(out) :: array(:, :, :, :, :)
425  character(len=*), intent(in) :: dsetname
426  integer, intent(out) :: error
427 
428 #define DATATYPE H5T_NATIVE_DOUBLE
429 #include "sll_k_hdf5_par_read_array.F90"
430 #undef DATATYPE
431 
432  end subroutine sll_hdf5_par_read_dble_array_5d
433 
434  !-----------------------------------------------------------------------------
436  !-----------------------------------------------------------------------------
438  handle, global_size, offset, array, dsetname, error)
439  integer, parameter :: dspace_dims = 6
440  type(sll_t_hdf5_par_handle), intent(in) :: handle
441  integer(i64), intent(in) :: global_size(dspace_dims)
442  integer(i64), intent(in) :: offset(dspace_dims)
443  real(f64), intent(out) :: array(:, :, :, :, :, :)
444  character(len=*), intent(in) :: dsetname
445  integer, intent(out) :: error
446 
447 #define DATATYPE H5T_NATIVE_DOUBLE
448 #include "sll_k_hdf5_par_read_array.F90"
449 #undef DATATYPE
450 
451  end subroutine sll_hdf5_par_read_dble_array_6d
452 
453 #endif
454 
455 end module sll_m_hdf5_io_parallel
Collectively read distributed nD array from HDF5 file.
Collectively write distributed nD array into HDF5 file.
Parallel version of sll_hdf5_io.
subroutine sll_hdf5_par_write_dble_array_3d(handle, global_size, offset, array, dsetname, error, chunk_dims)
Write 3D array of double precision floats into HDF5 file.
subroutine, public sll_s_hdf5_par_file_close(handle, error)
Close existing HDF5 file.
subroutine sll_hdf5_par_write_dble_array_1d(handle, global_size, offset, array, dsetname, error, chunk_dims)
Write 1D array of double precision floats into HDF5 file.
subroutine sll_hdf5_par_write_dble_array_2d(handle, global_size, offset, array, dsetname, error, chunk_dims)
Write 2D array of double precision floats into HDF5 file.
subroutine, public sll_s_hdf5_par_file_create(filename, comm, handle, error)
Create new HDF5 file.
subroutine sll_hdf5_par_write_dble_array_4d(handle, global_size, offset, array, dsetname, error, chunk_dims)
Write 4D array of double precision floats into HDF5 file.
subroutine, public sll_s_hdf5_par_file_open(filename, comm, handle, error)
Open existing HDF5 file.
subroutine sll_hdf5_par_read_dble_array_6d(handle, global_size, offset, array, dsetname, error)
Read 6D array of double precision floats from HDF5 file.
subroutine sll_hdf5_par_read_dble_array_5d(handle, global_size, offset, array, dsetname, error)
Read 5D array of double precision floats from HDF5 file.
subroutine sll_hdf5_par_write_dble_array_6d(handle, global_size, offset, array, dsetname, error, chunk_dims)
Write 6D array of double precision floats into HDF5 file.
subroutine sll_hdf5_par_read_dble_array_4d(handle, global_size, offset, array, dsetname, error)
Read 4D array of double precision floats from HDF5 file.
subroutine sll_hdf5_par_write_dble_array_5d(handle, global_size, offset, array, dsetname, error, chunk_dims)
Write 5D array of double precision floats into HDF5 file.
subroutine sll_hdf5_par_read_dble_array_3d(handle, global_size, offset, array, dsetname, error)
Read 3D array of double precision floats from HDF5 file.
subroutine sll_hdf5_par_read_dble_array_1d(handle, global_size, offset, array, dsetname, error)
Read 1D array of double precision floats from HDF5 file.
subroutine sll_hdf5_par_read_dble_array_2d(handle, global_size, offset, array, dsetname, error)
Read 2D array of double precision floats from HDF5 file.
    Report Typos and Errors