29 #include "sll_assert.h"
30 #include "sll_working_precision.h"
44 h5fd_mpio_collective_f, &
47 h5p_dataset_create_f, &
59 h5sselect_hyperslab_f, &
83 integer(hid_t),
private :: file_id
148 character(len=*),
intent(in) :: filename
149 integer,
intent(in) :: comm
151 integer,
intent(out) :: error
153 integer(hid_t) :: plist_id
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)
181 character(len=*),
intent(in) :: filename
182 integer,
intent(in) :: comm
184 integer,
intent(out) :: error
186 integer(hid_t) :: plist_id
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)
211 integer,
intent(out) :: error
213 call h5fclose_f(handle%file_id, error)
214 sll_assert_always(error == 0)
215 call h5close_f(error)
216 sll_assert_always(error == 0)
223 handle, global_size, offset, array, dsetname, error, chunk_dims)
224 integer,
parameter :: dspace_dims = 1
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)
233 #define DATATYPE H5T_NATIVE_DOUBLE
234 #include "sll_k_hdf5_par_write_array.F90"
243 handle, global_size, offset, array, dsetname, error, chunk_dims)
244 integer,
parameter :: dspace_dims = 2
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)
253 #define DATATYPE H5T_NATIVE_DOUBLE
254 #include "sll_k_hdf5_par_write_array.F90"
263 handle, global_size, offset, array, dsetname, error, chunk_dims)
264 integer,
parameter :: dspace_dims = 3
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)
273 #define DATATYPE H5T_NATIVE_DOUBLE
274 #include "sll_k_hdf5_par_write_array.F90"
283 handle, global_size, offset, array, dsetname, error, chunk_dims)
284 integer,
parameter :: dspace_dims = 4
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)
293 #define DATATYPE H5T_NATIVE_DOUBLE
294 #include "sll_k_hdf5_par_write_array.F90"
303 handle, global_size, offset, array, dsetname, error, chunk_dims)
304 integer,
parameter :: dspace_dims = 5
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)
313 #define DATATYPE H5T_NATIVE_DOUBLE
314 #include "sll_k_hdf5_par_write_array.F90"
323 handle, global_size, offset, array, dsetname, error, chunk_dims)
324 integer,
parameter :: dspace_dims = 6
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)
333 #define DATATYPE H5T_NATIVE_DOUBLE
334 #include "sll_k_hdf5_par_write_array.F90"
343 handle, global_size, offset, array, dsetname, error)
344 integer,
parameter :: dspace_dims = 1
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
352 #define DATATYPE H5T_NATIVE_DOUBLE
353 #include "sll_k_hdf5_par_read_array.F90"
362 handle, global_size, offset, array, dsetname, error)
363 integer,
parameter :: dspace_dims = 2
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
371 #define DATATYPE H5T_NATIVE_DOUBLE
372 #include "sll_k_hdf5_par_read_array.F90"
381 handle, global_size, offset, array, dsetname, error)
382 integer,
parameter :: dspace_dims = 3
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
390 #define DATATYPE H5T_NATIVE_DOUBLE
391 #include "sll_k_hdf5_par_read_array.F90"
400 handle, global_size, offset, array, dsetname, error)
401 integer,
parameter :: dspace_dims = 4
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
409 #define DATATYPE H5T_NATIVE_DOUBLE
410 #include "sll_k_hdf5_par_read_array.F90"
419 handle, global_size, offset, array, dsetname, error)
420 integer,
parameter :: dspace_dims = 5
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
428 #define DATATYPE H5T_NATIVE_DOUBLE
429 #include "sll_k_hdf5_par_read_array.F90"
438 handle, global_size, offset, array, dsetname, error)
439 integer,
parameter :: dspace_dims = 6
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
447 #define DATATYPE H5T_NATIVE_DOUBLE
448 #include "sll_k_hdf5_par_read_array.F90"
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.
Opaque object around HDF5 file id.