DAMASK with grid solvers  Revision: v2.0.3-2204-gdb1f2151
The Düsseldorf Advanced Material Simulation Kit with Grid Solvers
HDF5_utilities.f90
Go to the documentation of this file.
1 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/HDF5_utilities.f90"
2 # 1 "<built-in>"
3 # 1 "<command-line>"
4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/HDF5_utilities.f90"
5 !--------------------------------------------------------------------------------------------------
10 !--------------------------------------------------------------------------------------------------
12  use hdf5
13 
14  use petsc
15 
16 
17  use prec
18  use io
19  use rotations
20  use numerics
21 
22  implicit none
23  public
24 
25 !--------------------------------------------------------------------------------------------------
28 !--------------------------------------------------------------------------------------------------
29  interface hdf5_read
30  module procedure hdf5_read_real1
31  module procedure hdf5_read_real2
32  module procedure hdf5_read_real3
33  module procedure hdf5_read_real4
34  module procedure hdf5_read_real5
35  module procedure hdf5_read_real6
36  module procedure hdf5_read_real7
37 
38  module procedure hdf5_read_int1
39  module procedure hdf5_read_int2
40  module procedure hdf5_read_int3
41  module procedure hdf5_read_int4
42  module procedure hdf5_read_int5
43  module procedure hdf5_read_int6
44  module procedure hdf5_read_int7
45 
46  end interface hdf5_read
47 
48 !--------------------------------------------------------------------------------------------------
51 !--------------------------------------------------------------------------------------------------
52  interface hdf5_write
53  module procedure hdf5_write_real1
54  module procedure hdf5_write_real2
55  module procedure hdf5_write_real3
56  module procedure hdf5_write_real4
57  module procedure hdf5_write_real5
58  module procedure hdf5_write_real6
59  module procedure hdf5_write_real7
60 
61  module procedure hdf5_write_int1
62  module procedure hdf5_write_int2
63  module procedure hdf5_write_int3
64  module procedure hdf5_write_int4
65  module procedure hdf5_write_int5
66  module procedure hdf5_write_int6
67  module procedure hdf5_write_int7
68 
69  module procedure hdf5_write_rotation
70 
71  end interface hdf5_write
72 
73 !--------------------------------------------------------------------------------------------------
75 !--------------------------------------------------------------------------------------------------
77  module procedure hdf5_addattribute_str
78  module procedure hdf5_addattribute_int
79  module procedure hdf5_addattribute_real
80  module procedure hdf5_addattribute_int_array
81  module procedure hdf5_addattribute_real_array
82  end interface hdf5_addattribute
83 
84 contains
85 
86 
87 !--------------------------------------------------------------------------------------------------
89 !--------------------------------------------------------------------------------------------------
90 subroutine hdf5_utilities_init
91 
92  integer :: hdferr
93  integer(SIZE_T) :: typeSize
94 
95  write(6,'(/,a)') ' <<<+- HDF5_Utilities init -+>>>'
96 
97 !--------------------------------------------------------------------------------------------------
98 !initialize HDF5 library and check if integer and float type size match
99  call h5open_f(hdferr)
100  if (hdferr < 0) call io_error(1,ext_msg='HDF5_Utilities_init: h5open_f')
101 
102  call h5tget_size_f(h5t_native_integer,typesize, hdferr)
103  if (hdferr < 0) call io_error(1,ext_msg='HDF5_Utilities_init: h5tget_size_f (int)')
104  if (int(bit_size(0),size_t)/=typesize*8) &
105  call io_error(0,ext_msg='Default integer size does not match H5T_NATIVE_INTEGER')
106 
107  call h5tget_size_f(h5t_native_double,typesize, hdferr)
108  if (hdferr < 0) call io_error(1,ext_msg='HDF5_Utilities_init: h5tget_size_f (double)')
109  if (int(storage_size(0.0_preal),size_t)/=typesize*8) &
110  call io_error(0,ext_msg='pReal does not match H5T_NATIVE_DOUBLE')
111 
112 end subroutine hdf5_utilities_init
113 
114 
115 !--------------------------------------------------------------------------------------------------
117 !--------------------------------------------------------------------------------------------------
118 integer(HID_T) function hdf5_openfile(fileName,mode,parallel)
119 
120  character(len=*), intent(in) :: filename
121  character, intent(in), optional :: mode
122  logical, intent(in), optional :: parallel
123 
124  character :: m
125  integer(HID_T) :: plist_id
126  integer :: hdferr
127 
128  if (present(mode)) then
129  m = mode
130  else
131  m = 'r'
132  endif
133 
134  call h5pcreate_f(h5p_file_access_f, plist_id, hdferr)
135  if (hdferr < 0) call io_error(1,ext_msg='HDF5_openFile: h5pcreate_f')
136 
137 
138  if (present(parallel)) then; if (parallel) then
139  call h5pset_fapl_mpio_f(plist_id, petsc_comm_world, mpi_info_null, hdferr)
140  if (hdferr < 0) call io_error(1,ext_msg='HDF5_openFile: h5pset_fapl_mpio_f')
141  endif; endif
142 
143 
144  if (m == 'w') then
145  call h5fcreate_f(filename,h5f_acc_trunc_f,hdf5_openfile,hdferr,access_prp = plist_id)
146  if (hdferr < 0) call io_error(1,ext_msg='HDF5_openFile: h5fcreate_f (w)')
147  elseif(m == 'a') then
148  call h5fopen_f(filename,h5f_acc_rdwr_f,hdf5_openfile,hdferr,access_prp = plist_id)
149  if (hdferr < 0) call io_error(1,ext_msg='HDF5_openFile: h5fopen_f (a)')
150  elseif(m == 'r') then
151  call h5fopen_f(filename,h5f_acc_rdonly_f,hdf5_openfile,hdferr,access_prp = plist_id)
152  if (hdferr < 0) call io_error(1,ext_msg='HDF5_openFile: h5fopen_f (r)')
153  else
154  call io_error(1,ext_msg='HDF5_openFile: h5fopen_f unknown access mode: '//trim(m))
155  endif
156 
157  call h5pclose_f(plist_id, hdferr)
158  if (hdferr < 0) call io_error(1,ext_msg='HDF5_openFile: h5pclose_f')
159 
160 end function hdf5_openfile
161 
162 
163 !--------------------------------------------------------------------------------------------------
165 !--------------------------------------------------------------------------------------------------
166 subroutine hdf5_closefile(fileHandle)
167 
168  integer(HID_T), intent(in) :: fileHandle
169 
170  integer :: hdferr
171 
172  call h5fclose_f(filehandle,hdferr)
173  if (hdferr < 0) call io_error(1,ext_msg='HDF5_closeFile: h5fclose_f')
174 
175 end subroutine hdf5_closefile
176 
177 
178 !--------------------------------------------------------------------------------------------------
180 !--------------------------------------------------------------------------------------------------
181 integer(HID_T) function hdf5_addgroup(fileHandle,groupName)
182 
183  integer(HID_T), intent(in) :: filehandle
184  character(len=*), intent(in) :: groupname
185 
186  integer :: hdferr
187  integer(HID_T) :: aplist_id
188 
189 !-------------------------------------------------------------------------------------------------
190 ! creating a property list for data access properties
191  call h5pcreate_f(h5p_group_access_f, aplist_id, hdferr)
192  if (hdferr < 0) call io_error(1,ext_msg = 'HDF5_addGroup: h5pcreate_f ('//trim(groupname)//')')
193 
194 !-------------------------------------------------------------------------------------------------
195 ! setting I/O mode to collective
196 
197  call h5pset_all_coll_metadata_ops_f(aplist_id, .true., hdferr)
198  if (hdferr < 0) call io_error(1,ext_msg = 'HDF5_addGroup: h5pset_all_coll_metadata_ops_f ('//trim(groupname)//')')
199 
200 
201 !-------------------------------------------------------------------------------------------------
202 ! Create group
203  call h5gcreate_f(filehandle, trim(groupname), hdf5_addgroup, hdferr, object_namelen_default_f,gapl_id = aplist_id)
204  if (hdferr < 0) call io_error(1,ext_msg = 'HDF5_addGroup: h5gcreate_f ('//trim(groupname)//')')
205 
206  call h5pclose_f(aplist_id,hdferr)
207 
208 end function hdf5_addgroup
209 
210 
211 !--------------------------------------------------------------------------------------------------
213 !--------------------------------------------------------------------------------------------------
214 integer(HID_T) function hdf5_opengroup(fileHandle,groupName)
215 
216  integer(HID_T), intent(in) :: filehandle
217  character(len=*), intent(in) :: groupname
218 
219 
220  integer :: hdferr
221  integer(HID_T) :: aplist_id
222  logical :: is_collective
223 
224 
225  !-------------------------------------------------------------------------------------------------
226  ! creating a property list for data access properties
227  call h5pcreate_f(h5p_group_access_f, aplist_id, hdferr)
228  if (hdferr < 0) call io_error(1,ext_msg = 'HDF5_openGroup: h5pcreate_f ('//trim(groupname)//')')
229 
230  !-------------------------------------------------------------------------------------------------
231  ! setting I/O mode to collective
232 
233  call h5pget_all_coll_metadata_ops_f(aplist_id, is_collective, hdferr)
234  if (hdferr < 0) call io_error(1,ext_msg = 'HDF5_openGroup: h5pset_all_coll_metadata_ops_f ('//trim(groupname)//')')
235 
236 
237  !-------------------------------------------------------------------------------------------------
238  ! opening the group
239  call h5gopen_f(filehandle, trim(groupname), hdf5_opengroup, hdferr, gapl_id = aplist_id)
240  if (hdferr < 0) call io_error(1,ext_msg = 'HDF5_openGroup: h5gopen_f ('//trim(groupname)//')')
241 
242  call h5pclose_f(aplist_id,hdferr)
243 
244 end function hdf5_opengroup
245 
246 
247 !--------------------------------------------------------------------------------------------------
249 !--------------------------------------------------------------------------------------------------
250 subroutine hdf5_closegroup(group_id)
251 
252  integer(HID_T), intent(in) :: group_id
253  integer :: hdferr
254 
255  call h5gclose_f(group_id, hdferr)
256  if (hdferr < 0) call io_error(1,ext_msg = 'HDF5_closeGroup: h5gclose_f (el is ID)', el = int(group_id))
257 
258 end subroutine hdf5_closegroup
259 
260 
261 !--------------------------------------------------------------------------------------------------
263 !--------------------------------------------------------------------------------------------------
264 logical function hdf5_objectexists(loc_id,path)
265 
266  integer(HID_T), intent(in) :: loc_id
267  character(len=*), intent(in), optional :: path
268 
269  integer :: hdferr
270  character(len=pStringLen) :: p
271 
272  if (present(path)) then
273  p = trim(path)
274  else
275  p = '.'
276  endif
277 
278  call h5lexists_f(loc_id, p, hdf5_objectexists, hdferr)
279  if (hdferr < 0) call io_error(1,ext_msg = 'HDF5_objectExists: h5oexists_by_name_f')
280 
281  if(hdf5_objectexists) then
282  call h5oexists_by_name_f(loc_id, p, hdf5_objectexists, hdferr)
283  if (hdferr < 0) call io_error(1,ext_msg = 'HDF5_objectExists: h5oexists_by_name_f')
284  endif
285 
286 end function hdf5_objectexists
287 
288 
289 !--------------------------------------------------------------------------------------------------
291 !--------------------------------------------------------------------------------------------------
292 subroutine hdf5_addattribute_str(loc_id,attrLabel,attrValue,path)
293 
294  integer(HID_T), intent(in) :: loc_id
295  character(len=*), intent(in) :: attrLabel, attrValue
296  character(len=*), intent(in), optional :: path
297 
298  integer :: hdferr
299  integer(HID_T) :: attr_id, space_id, type_id
300  logical :: attrExists
301  character(len=pStringLen) :: p
302 
303  if (present(path)) then
304  p = trim(path)
305  else
306  p = '.'
307  endif
308 
309  call h5screate_f(h5s_scalar_f,space_id,hdferr)
310  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_str: h5screate_f')
311  call h5tcopy_f(h5t_native_character, type_id, hdferr)
312  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_str: h5tcopy_f')
313  call h5tset_size_f(type_id, int(len_trim(attrvalue),hsize_t), hdferr)
314  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_str: h5tset_size_f')
315  call h5aexists_by_name_f(loc_id,trim(p),attrlabel,attrexists,hdferr)
316  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_str: h5aexists_by_name_f')
317  if (attrexists) then
318  call h5adelete_by_name_f(loc_id, trim(p), attrlabel, hdferr)
319  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_str: h5adelete_by_name_f')
320  endif
321  call h5acreate_by_name_f(loc_id,trim(p),trim(attrlabel),type_id,space_id,attr_id,hdferr)
322  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_str: h5acreate_f')
323  call h5awrite_f(attr_id, type_id, trim(attrvalue), int([1],hsize_t), hdferr)
324  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_str: h5awrite_f')
325  call h5aclose_f(attr_id,hdferr)
326  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_str: h5aclose_f')
327  call h5tclose_f(type_id,hdferr)
328  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_str: h5tclose_f')
329  call h5sclose_f(space_id,hdferr)
330  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_str: h5sclose_f')
331 
332 end subroutine hdf5_addattribute_str
333 
334 
335 !--------------------------------------------------------------------------------------------------
337 !--------------------------------------------------------------------------------------------------
338 subroutine hdf5_addattribute_int(loc_id,attrLabel,attrValue,path)
339 
340  integer(HID_T), intent(in) :: loc_id
341  character(len=*), intent(in) :: attrLabel
342  integer, intent(in) :: attrValue
343  character(len=*), intent(in), optional :: path
344 
345  integer :: hdferr
346  integer(HID_T) :: attr_id, space_id
347  logical :: attrExists
348  character(len=pStringLen) :: p
349 
350  if (present(path)) then
351  p = trim(path)
352  else
353  p = '.'
354  endif
355 
356  call h5screate_f(h5s_scalar_f,space_id,hdferr)
357  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int: h5screate_f')
358  call h5aexists_by_name_f(loc_id,trim(p),attrlabel,attrexists,hdferr)
359  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int: h5aexists_by_name_f')
360  if (attrexists) then
361  call h5adelete_by_name_f(loc_id, trim(p), attrlabel, hdferr)
362  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int: h5adelete_by_name_f')
363  endif
364  call h5acreate_by_name_f(loc_id,trim(p),trim(attrlabel),h5t_native_integer,space_id,attr_id,hdferr)
365  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int: h5acreate_f')
366  call h5awrite_f(attr_id, h5t_native_integer, attrvalue, int([1],hsize_t), hdferr)
367  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int: h5awrite_f')
368  call h5aclose_f(attr_id,hdferr)
369  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int: h5tclose_f')
370  call h5sclose_f(space_id,hdferr)
371  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int: h5sclose_f')
372 
373 end subroutine hdf5_addattribute_int
374 
375 
376 !--------------------------------------------------------------------------------------------------
378 !--------------------------------------------------------------------------------------------------
379 subroutine hdf5_addattribute_real(loc_id,attrLabel,attrValue,path)
380 
381  integer(HID_T), intent(in) :: loc_id
382  character(len=*), intent(in) :: attrLabel
383  real(pReal), intent(in) :: attrValue
384  character(len=*), intent(in), optional :: path
385 
386  integer :: hdferr
387  integer(HID_T) :: attr_id, space_id
388  logical :: attrExists
389  character(len=pStringLen) :: p
390 
391  if (present(path)) then
392  p = trim(path)
393  else
394  p = '.'
395  endif
396 
397  call h5screate_f(h5s_scalar_f,space_id,hdferr)
398  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_real: h5screate_f')
399  call h5aexists_by_name_f(loc_id,trim(p),attrlabel,attrexists,hdferr)
400  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_real: h5aexists_by_name_f')
401  if (attrexists) then
402  call h5adelete_by_name_f(loc_id, trim(p), attrlabel, hdferr)
403  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_real: h5adelete_by_name_f')
404  endif
405  call h5acreate_by_name_f(loc_id,trim(p),trim(attrlabel),h5t_native_double,space_id,attr_id,hdferr)
406  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_real: h5acreate_f')
407  call h5awrite_f(attr_id, h5t_native_double, attrvalue, int([1],hsize_t), hdferr)
408  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_real: h5awrite_f')
409  call h5aclose_f(attr_id,hdferr)
410  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_real: h5tclose_f')
411  call h5sclose_f(space_id,hdferr)
412  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_real: h5sclose_f')
413 
414 end subroutine hdf5_addattribute_real
415 
416 
417 !--------------------------------------------------------------------------------------------------
419 !--------------------------------------------------------------------------------------------------
420 subroutine hdf5_addattribute_int_array(loc_id,attrLabel,attrValue,path)
421 
422  integer(HID_T), intent(in) :: loc_id
423  character(len=*), intent(in) :: attrLabel
424  integer, intent(in), dimension(:) :: attrValue
425  character(len=*), intent(in), optional :: path
426 
427  integer :: hdferr
428  integer(HID_T) :: attr_id, space_id
429  integer(HSIZE_T),dimension(1) :: array_size
430  logical :: attrExists
431  character(len=pStringLen) :: p
432 
433  if (present(path)) then
434  p = trim(path)
435  else
436  p = '.'
437  endif
438 
439  array_size = size(attrvalue,kind=hsize_t)
440 
441  call h5screate_simple_f(1, array_size, space_id, hdferr, array_size)
442  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int_array: h5screate_f')
443  call h5aexists_by_name_f(loc_id,trim(p),attrlabel,attrexists,hdferr)
444  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int_array: h5aexists_by_name_f')
445  if (attrexists) then
446  call h5adelete_by_name_f(loc_id, trim(p), attrlabel, hdferr)
447  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int_array: h5adelete_by_name_f')
448  endif
449  call h5acreate_by_name_f(loc_id,trim(p),trim(attrlabel),h5t_native_integer,space_id,attr_id,hdferr)
450  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int_array: h5acreate_f')
451  call h5awrite_f(attr_id, h5t_native_integer, attrvalue, array_size, hdferr)
452  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int_array: h5awrite_f')
453  call h5aclose_f(attr_id,hdferr)
454  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int_array: h5tclose_f')
455  call h5sclose_f(space_id,hdferr)
456  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int_array: h5sclose_f')
457 
458 end subroutine hdf5_addattribute_int_array
459 
460 
461 !--------------------------------------------------------------------------------------------------
463 !--------------------------------------------------------------------------------------------------
464 subroutine hdf5_addattribute_real_array(loc_id,attrLabel,attrValue,path)
465 
466  integer(HID_T), intent(in) :: loc_id
467  character(len=*), intent(in) :: attrLabel
468  real(pReal), intent(in), dimension(:) :: attrValue
469  character(len=*), intent(in), optional :: path
470 
471  integer :: hdferr
472  integer(HID_T) :: attr_id, space_id
473  integer(HSIZE_T),dimension(1) :: array_size
474  logical :: attrExists
475  character(len=pStringLen) :: p
476 
477  if (present(path)) then
478  p = trim(path)
479  else
480  p = '.'
481  endif
482 
483  array_size = size(attrvalue,kind=hsize_t)
484 
485  call h5screate_simple_f(1, array_size, space_id, hdferr, array_size)
486  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int_array: h5screate_f')
487  call h5aexists_by_name_f(loc_id,trim(p),attrlabel,attrexists,hdferr)
488  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int_array: h5aexists_by_name_f')
489  if (attrexists) then
490  call h5adelete_by_name_f(loc_id, trim(p), attrlabel, hdferr)
491  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int_array: h5adelete_by_name_f')
492  endif
493  call h5acreate_by_name_f(loc_id,trim(p),trim(attrlabel),h5t_native_double,space_id,attr_id,hdferr)
494  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int_array: h5acreate_f')
495  call h5awrite_f(attr_id, h5t_native_double, attrvalue, array_size, hdferr)
496  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int_array: h5awrite_f')
497  call h5aclose_f(attr_id,hdferr)
498  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int_array: h5tclose_f')
499  call h5sclose_f(space_id,hdferr)
500  if (hdferr < 0) call io_error(1,ext_msg='HDF5_addAttribute_int_array: h5sclose_f')
501 
502 end subroutine hdf5_addattribute_real_array
503 
504 
505 !--------------------------------------------------------------------------------------------------
507 !--------------------------------------------------------------------------------------------------
508 subroutine hdf5_setlink(loc_id,target_name,link_name)
509 
510  character(len=*), intent(in) :: target_name, link_name
511  integer(HID_T), intent(in) :: loc_id
512  integer :: hdferr
513  logical :: linkExists
514 
515  call h5lexists_f(loc_id, link_name,linkexists, hdferr)
516  if (hdferr < 0) call io_error(1,ext_msg = 'HDF5_setLink: h5lexists_soft_f ('//trim(link_name)//')')
517  if (linkexists) then
518  call h5ldelete_f(loc_id,link_name, hdferr)
519  if (hdferr < 0) call io_error(1,ext_msg = 'HDF5_setLink: h5ldelete_soft_f ('//trim(link_name)//')')
520  endif
521  call h5lcreate_soft_f(target_name, loc_id, link_name, hdferr)
522  if (hdferr < 0) call io_error(1,ext_msg = 'HDF5_setLink: h5lcreate_soft_f ('//trim(target_name)//' '//trim(link_name)//')')
523 
524 end subroutine hdf5_setlink
525 
526 
527 !--------------------------------------------------------------------------------------------------
529 !--------------------------------------------------------------------------------------------------
530 subroutine hdf5_read_real1(loc_id,dataset,datasetName,parallel)
531 
532  real(pReal), intent(out), dimension(:) :: dataset
533  integer(HID_T), intent(in) :: loc_id
534  character(len=*), intent(in) :: datasetName
535  logical, intent(in), optional :: parallel
536 
537  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id, aplist_id
538  integer(HSIZE_T), dimension(size(shape(dataset))) :: & ! ToDo: Fortran 2018 size(shape(A)) = rank(A)
539  myStart, &
540  myShape, & !< shape of the dataset (this process)
541  totalShape
542  integer :: hdferr
543 
544 !---------------------------------------------------------------------------------------------------
545 ! determine shape of dataset
546  myshape = int(shape(dataset),hsize_t)
547  if (any(myshape(1:size(myshape)-1) == 0)) return
548 
549 !---------------------------------------------------------------------------------------------------
550 ! initialize HDF5 data structures
551  if (present(parallel)) then
552  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
553  mystart, totalshape, loc_id,myshape,datasetname,parallel)
554  else
555  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
556  mystart, totalshape, loc_id,myshape,datasetname,.false.)
557  endif
558 
559  call h5dread_f(dset_id, h5t_native_double,dataset,totalshape, hdferr,&
560  file_space_id = filespace_id, xfer_prp = plist_id, mem_space_id = memspace_id)
561  if (hdferr < 0) call io_error(1,ext_msg='HDF5_read_real1: h5dread_f')
562 
563  call finalize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id)
564 
565 end subroutine hdf5_read_real1
566 
567 !--------------------------------------------------------------------------------------------------
569 !--------------------------------------------------------------------------------------------------
570 subroutine hdf5_read_real2(loc_id,dataset,datasetName,parallel)
571 
572  real(pReal), intent(out), dimension(:,:) :: dataset
573  integer(HID_T), intent(in) :: loc_id
574  character(len=*), intent(in) :: datasetName
575  logical, intent(in), optional :: parallel
576 
577  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id, aplist_id
578  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
579  myStart, &
580  myShape, & !< shape of the dataset (this process)
581  totalShape
582  integer :: hdferr
583 
584 !---------------------------------------------------------------------------------------------------
585 ! determine shape of dataset
586  myshape = int(shape(dataset),hsize_t)
587  if (any(myshape(1:size(myshape)-1) == 0)) return
588 
589 !---------------------------------------------------------------------------------------------------
590 ! initialize HDF5 data structures
591  if (present(parallel)) then
592  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
593  mystart, totalshape, loc_id,myshape,datasetname,parallel)
594  else
595  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
596  mystart, totalshape, loc_id,myshape,datasetname,.false.)
597  endif
598 
599  call h5dread_f(dset_id, h5t_native_double,dataset,totalshape, hdferr,&
600  file_space_id = filespace_id, xfer_prp = plist_id, mem_space_id = memspace_id)
601  if (hdferr < 0) call io_error(1,ext_msg='HDF5_read_real2: h5dread_f')
602 
603  call finalize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id)
604 
605 end subroutine hdf5_read_real2
606 
607 !--------------------------------------------------------------------------------------------------
609 !--------------------------------------------------------------------------------------------------
610 subroutine hdf5_read_real3(loc_id,dataset,datasetName,parallel)
611 
612  real(pReal), intent(out), dimension(:,:,:) :: dataset
613  integer(HID_T), intent(in) :: loc_id
614  character(len=*), intent(in) :: datasetName
615  logical, intent(in), optional :: parallel
616 
617  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id, aplist_id
618  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
619  myStart, &
620  myShape, & !< shape of the dataset (this process)
621  totalShape
622  integer :: hdferr
623 
624 !---------------------------------------------------------------------------------------------------
625 ! determine shape of dataset
626  myshape = int(shape(dataset),hsize_t)
627  if (any(myshape(1:size(myshape)-1) == 0)) return
628 
629 !---------------------------------------------------------------------------------------------------
630 ! initialize HDF5 data structures
631  if (present(parallel)) then
632  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
633  mystart, totalshape, loc_id,myshape,datasetname,parallel)
634  else
635  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
636  mystart, totalshape, loc_id,myshape,datasetname,.false.)
637  endif
638 
639  call h5dread_f(dset_id, h5t_native_double,dataset,totalshape, hdferr,&
640  file_space_id = filespace_id, xfer_prp = plist_id, mem_space_id = memspace_id)
641  if (hdferr < 0) call io_error(1,ext_msg='HDF5_read_real3: h5dread_f')
642 
643  call finalize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id)
644 
645 end subroutine hdf5_read_real3
646 
647 !--------------------------------------------------------------------------------------------------
649 !--------------------------------------------------------------------------------------------------
650 subroutine hdf5_read_real4(loc_id,dataset,datasetName,parallel)
651 
652  real(pReal), intent(out), dimension(:,:,:,:) :: dataset
653  integer(HID_T), intent(in) :: loc_id
654  character(len=*), intent(in) :: datasetName
655  logical, intent(in), optional :: parallel
656 
657  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id, aplist_id
658  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
659  myStart, &
660  myShape, & !< shape of the dataset (this process)
661  totalShape
662  integer :: hdferr
663 
664 !---------------------------------------------------------------------------------------------------
665 ! determine shape of dataset
666  myshape = int(shape(dataset),hsize_t)
667  if (any(myshape(1:size(myshape)-1) == 0)) return
668 
669 !---------------------------------------------------------------------------------------------------
670 ! initialize HDF5 data structures
671  if (present(parallel)) then
672  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
673  mystart, totalshape, loc_id,myshape,datasetname,parallel)
674  else
675  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
676  mystart, totalshape, loc_id,myshape,datasetname,.false.)
677  endif
678 
679  call h5dread_f(dset_id, h5t_native_double,dataset,totalshape, hdferr,&
680  file_space_id = filespace_id, xfer_prp = plist_id, mem_space_id = memspace_id)
681  if (hdferr < 0) call io_error(1,ext_msg='HDF5_read_real4: h5dread_f')
682 
683  call finalize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id)
684 
685 end subroutine hdf5_read_real4
686 
687 !--------------------------------------------------------------------------------------------------
689 !--------------------------------------------------------------------------------------------------
690 subroutine hdf5_read_real5(loc_id,dataset,datasetName,parallel)
691 
692  real(pReal), intent(out), dimension(:,:,:,:,:) :: dataset
693  integer(HID_T), intent(in) :: loc_id
694  character(len=*), intent(in) :: datasetName
695  logical, intent(in), optional :: parallel
696 
697  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id, aplist_id
698  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
699  myStart, &
700  myShape, & !< shape of the dataset (this process)
701  totalShape
702  integer :: hdferr
703 
704 !---------------------------------------------------------------------------------------------------
705 ! determine shape of dataset
706  myshape = int(shape(dataset),hsize_t)
707  if (any(myshape(1:size(myshape)-1) == 0)) return
708 
709 !---------------------------------------------------------------------------------------------------
710 ! initialize HDF5 data structures
711  if (present(parallel)) then
712  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
713  mystart, totalshape, loc_id,myshape,datasetname,parallel)
714  else
715  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
716  mystart, totalshape, loc_id,myshape,datasetname,.false.)
717  endif
718 
719  call h5dread_f(dset_id, h5t_native_double,dataset,totalshape, hdferr,&
720  file_space_id = filespace_id, xfer_prp = plist_id, mem_space_id = memspace_id)
721  if (hdferr < 0) call io_error(1,ext_msg='HDF5_read_real5: h5dread_f')
722 
723  call finalize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id)
724 
725 end subroutine hdf5_read_real5
726 
727 !--------------------------------------------------------------------------------------------------
729 !--------------------------------------------------------------------------------------------------
730 subroutine hdf5_read_real6(loc_id,dataset,datasetName,parallel)
731 
732  real(pReal), intent(out), dimension(:,:,:,:,:,:) :: dataset
733  integer(HID_T), intent(in) :: loc_id
734  character(len=*), intent(in) :: datasetName
735  logical, intent(in), optional :: parallel
736 
737  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id, aplist_id
738  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
739  myStart, &
740  myShape, & !< shape of the dataset (this process)
741  totalShape
742  integer :: hdferr
743 
744 !---------------------------------------------------------------------------------------------------
745 ! determine shape of dataset
746  myshape = int(shape(dataset),hsize_t)
747  if (any(myshape(1:size(myshape)-1) == 0)) return
748 
749 !---------------------------------------------------------------------------------------------------
750 ! initialize HDF5 data structures
751  if (present(parallel)) then
752  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
753  mystart, totalshape, loc_id,myshape,datasetname,parallel)
754  else
755  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
756  mystart, totalshape, loc_id,myshape,datasetname,.false.)
757  endif
758 
759  call h5dread_f(dset_id, h5t_native_double,dataset,totalshape, hdferr,&
760  file_space_id = filespace_id, xfer_prp = plist_id, mem_space_id = memspace_id)
761  if (hdferr < 0) call io_error(1,ext_msg='HDF5_read_real6: h5dread_f')
762 
763  call finalize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id)
764 
765 end subroutine hdf5_read_real6
766 
767 !--------------------------------------------------------------------------------------------------
769 !--------------------------------------------------------------------------------------------------
770 subroutine hdf5_read_real7(loc_id,dataset,datasetName,parallel)
771 
772  real(pReal), intent(out), dimension(:,:,:,:,:,:,:) :: dataset
773  integer(HID_T), intent(in) :: loc_id
774  character(len=*), intent(in) :: datasetName
775  logical, intent(in), optional :: parallel
776 
777  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id, aplist_id
778  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
779  myStart, &
780  myShape, & !< shape of the dataset (this process)
781  totalShape
782  integer :: hdferr
783 
784 !---------------------------------------------------------------------------------------------------
785 ! determine shape of dataset
786  myshape = int(shape(dataset),hsize_t)
787  if (any(myshape(1:size(myshape)-1) == 0)) return
788 
789 !---------------------------------------------------------------------------------------------------
790 ! initialize HDF5 data structures
791  if (present(parallel)) then
792  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
793  mystart, totalshape, loc_id,myshape,datasetname,parallel)
794  else
795  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
796  mystart, totalshape, loc_id,myshape,datasetname,.false.)
797  endif
798 
799  call h5dread_f(dset_id, h5t_native_double,dataset,totalshape, hdferr,&
800  file_space_id = filespace_id, xfer_prp = plist_id, mem_space_id = memspace_id)
801  if (hdferr < 0) call io_error(1,ext_msg='HDF5_read_real7: h5dread_f')
802 
803  call finalize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id)
804 
805 end subroutine hdf5_read_real7
806 
807 
808 !--------------------------------------------------------------------------------------------------
810 !--------------------------------------------------------------------------------------------------
811 subroutine hdf5_read_int1(loc_id,dataset,datasetName,parallel)
812 
813  integer, intent(out), dimension(:) :: dataset
814  integer(HID_T), intent(in) :: loc_id
815  character(len=*), intent(in) :: datasetName
816  logical, intent(in), optional :: parallel
817 
818 
819  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id, aplist_id
820  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
821  myStart, &
822  myShape, & !< shape of the dataset (this process)
823  totalShape
824  integer :: hdferr
825 
826 !---------------------------------------------------------------------------------------------------
827 ! determine shape of dataset
828  myshape = int(shape(dataset),hsize_t)
829  if (any(myshape(1:size(myshape)-1) == 0)) return
830 
831 !---------------------------------------------------------------------------------------------------
832 ! initialize HDF5 data structures
833  if (present(parallel)) then
834  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
835  mystart, totalshape, loc_id,myshape,datasetname,parallel)
836  else
837  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
838  mystart, totalshape, loc_id,myshape,datasetname,.false.)
839  endif
840 
841  call h5dread_f(dset_id, h5t_native_integer,dataset,totalshape, hdferr,&
842  file_space_id = filespace_id, xfer_prp = plist_id, mem_space_id = memspace_id)
843  if (hdferr < 0) call io_error(1,ext_msg='HDF5_read_int1: h5dread_f')
844 
845  call finalize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id)
846 
847 end subroutine hdf5_read_int1
848 
849 !--------------------------------------------------------------------------------------------------
851 !--------------------------------------------------------------------------------------------------
852 subroutine hdf5_read_int2(loc_id,dataset,datasetName,parallel)
853 
854  integer, intent(out), dimension(:,:) :: dataset
855  integer(HID_T), intent(in) :: loc_id
856  character(len=*), intent(in) :: datasetName
857  logical, intent(in), optional :: parallel
858 
859  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id, aplist_id
860  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
861  myStart, &
862  myShape, & !< shape of the dataset (this process)
863  totalShape
864  integer :: hdferr
865 
866 !---------------------------------------------------------------------------------------------------
867 ! determine shape of dataset
868  myshape = int(shape(dataset),hsize_t)
869  if (any(myshape(1:size(myshape)-1) == 0)) return
870 
871 !---------------------------------------------------------------------------------------------------
872 ! initialize HDF5 data structures
873  if (present(parallel)) then
874  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
875  mystart, totalshape, loc_id,myshape,datasetname,parallel)
876  else
877  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
878  mystart, totalshape, loc_id,myshape,datasetname,.false.)
879  endif
880 
881  call h5dread_f(dset_id, h5t_native_integer,dataset,totalshape, hdferr,&
882  file_space_id = filespace_id, xfer_prp = plist_id, mem_space_id = memspace_id)
883  if (hdferr < 0) call io_error(1,ext_msg='HDF5_read_int2: h5dread_f')
884 
885  call finalize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id)
886 
887 end subroutine hdf5_read_int2
888 
889 !--------------------------------------------------------------------------------------------------
891 !--------------------------------------------------------------------------------------------------
892 subroutine hdf5_read_int3(loc_id,dataset,datasetName,parallel)
893 
894  integer, intent(out), dimension(:,:,:) :: dataset
895  integer(HID_T), intent(in) :: loc_id
896  character(len=*), intent(in) :: datasetName
897  logical, intent(in), optional :: parallel
898 
899  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id, aplist_id
900  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
901  myStart, &
902  myShape, & !< shape of the dataset (this process)
903  totalShape
904  integer :: hdferr
905 
906 !---------------------------------------------------------------------------------------------------
907 ! determine shape of dataset
908  myshape = int(shape(dataset),hsize_t)
909  if (any(myshape(1:size(myshape)-1) == 0)) return
910 
911 !---------------------------------------------------------------------------------------------------
912 ! initialize HDF5 data structures
913  if (present(parallel)) then
914  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
915  mystart, totalshape, loc_id,myshape,datasetname,parallel)
916  else
917  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
918  mystart, totalshape, loc_id,myshape,datasetname,.false.)
919  endif
920 
921  call h5dread_f(dset_id, h5t_native_integer,dataset,totalshape, hdferr,&
922  file_space_id = filespace_id, xfer_prp = plist_id, mem_space_id = memspace_id)
923  if (hdferr < 0) call io_error(1,ext_msg='HDF5_read_int3: h5dread_f')
924 
925  call finalize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id)
926 
927 end subroutine hdf5_read_int3
928 
929 !--------------------------------------------------------------------------------------------------
931 !--------------------------------------------------------------------------------------------------
932 subroutine hdf5_read_int4(loc_id,dataset,datasetName,parallel)
933 
934  integer, intent(out), dimension(:,:,:,:) :: dataset
935  integer(HID_T), intent(in) :: loc_id
936  character(len=*), intent(in) :: datasetName
937  logical, intent(in), optional :: parallel
938 
939  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id, aplist_id
940  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
941  myStart, &
942  myShape, & !< shape of the dataset (this process)
943  totalShape
944  integer :: hdferr
945 
946 !---------------------------------------------------------------------------------------------------
947 ! determine shape of dataset
948  myshape = int(shape(dataset),hsize_t)
949  if (any(myshape(1:size(myshape)-1) == 0)) return
950 
951 !---------------------------------------------------------------------------------------------------
952 ! initialize HDF5 data structures
953  if (present(parallel)) then
954  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
955  mystart, totalshape, loc_id,myshape,datasetname,parallel)
956  else
957  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
958  mystart, totalshape, loc_id,myshape,datasetname,.false.)
959  endif
960 
961  call h5dread_f(dset_id, h5t_native_integer,dataset,totalshape, hdferr,&
962  file_space_id = filespace_id, xfer_prp = plist_id, mem_space_id = memspace_id)
963  if (hdferr < 0) call io_error(1,ext_msg='HDF5_read_int4: h5dread_f')
964 
965  call finalize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id)
966 
967 end subroutine hdf5_read_int4
968 
969 !--------------------------------------------------------------------------------------------------
971 !--------------------------------------------------------------------------------------------------
972 subroutine hdf5_read_int5(loc_id,dataset,datasetName,parallel)
973 
974  integer, intent(out), dimension(:,:,:,:,:) :: dataset
975  integer(HID_T), intent(in) :: loc_id
976  character(len=*), intent(in) :: datasetName
977  logical, intent(in), optional :: parallel
978 
979  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id, aplist_id
980  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
981  myStart, &
982  myShape, & !< shape of the dataset (this process)
983  totalShape
984  integer :: hdferr
985 
986 !---------------------------------------------------------------------------------------------------
987 ! determine shape of dataset
988  myshape = int(shape(dataset),hsize_t)
989  if (any(myshape(1:size(myshape)-1) == 0)) return
990 
991 !---------------------------------------------------------------------------------------------------
992 ! initialize HDF5 data structures
993  if (present(parallel)) then
994  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
995  mystart, totalshape, loc_id,myshape,datasetname,parallel)
996  else
997  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
998  mystart, totalshape, loc_id,myshape,datasetname,.false.)
999  endif
1000 
1001  call h5dread_f(dset_id, h5t_native_integer,dataset,totalshape, hdferr,&
1002  file_space_id = filespace_id, xfer_prp = plist_id, mem_space_id = memspace_id)
1003  if (hdferr < 0) call io_error(1,ext_msg='HDF5_read_int5: h5dread_f')
1004 
1005  call finalize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id)
1006 
1007 end subroutine hdf5_read_int5
1008 
1009 !--------------------------------------------------------------------------------------------------
1011 !--------------------------------------------------------------------------------------------------
1012 subroutine hdf5_read_int6(loc_id,dataset,datasetName,parallel)
1014  integer, intent(out), dimension(:,:,:,:,:,:) :: dataset
1015  integer(HID_T), intent(in) :: loc_id
1016  character(len=*), intent(in) :: datasetName
1017  logical, intent(in), optional :: parallel
1018 
1019  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id, aplist_id
1020  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
1021  myStart, &
1022  myShape, & !< shape of the dataset (this process)
1023  totalShape
1024  integer :: hdferr
1025 
1026 !---------------------------------------------------------------------------------------------------
1027 ! determine shape of dataset
1028  myshape = int(shape(dataset),hsize_t)
1029  if (any(myshape(1:size(myshape)-1) == 0)) return
1030 
1031 !---------------------------------------------------------------------------------------------------
1032 ! initialize HDF5 data structures
1033  if (present(parallel)) then
1034  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
1035  mystart, totalshape, loc_id,myshape,datasetname,parallel)
1036  else
1037  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
1038  mystart, totalshape, loc_id,myshape,datasetname,.false.)
1039  endif
1040 
1041  call h5dread_f(dset_id, h5t_native_integer,dataset,totalshape, hdferr,&
1042  file_space_id = filespace_id, xfer_prp = plist_id, mem_space_id = memspace_id)
1043  if (hdferr < 0) call io_error(1,ext_msg='HDF5_read_int6: h5dread_f')
1044 
1045  call finalize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id)
1046 
1047 end subroutine hdf5_read_int6
1048 
1049 !--------------------------------------------------------------------------------------------------
1051 !--------------------------------------------------------------------------------------------------
1052 subroutine hdf5_read_int7(loc_id,dataset,datasetName,parallel)
1054  integer, intent(out), dimension(:,:,:,:,:,:,:) :: dataset
1055  integer(HID_T), intent(in) :: loc_id
1056  character(len=*), intent(in) :: datasetName
1057  logical, intent(in), optional :: parallel
1058 
1059  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id, aplist_id
1060  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
1061  myStart, &
1062  myShape, & !< shape of the dataset (this process)
1063  totalShape
1064  integer :: hdferr
1065 
1066 !---------------------------------------------------------------------------------------------------
1067 ! determine shape of dataset
1068  myshape = int(shape(dataset),hsize_t)
1069  if (any(myshape(1:size(myshape)-1) == 0)) return
1070 
1071 !---------------------------------------------------------------------------------------------------
1072 ! initialize HDF5 data structures
1073  if (present(parallel)) then
1074  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
1075  mystart, totalshape, loc_id,myshape,datasetname,parallel)
1076  else
1077  call initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
1078  mystart, totalshape, loc_id,myshape,datasetname,.false.)
1079  endif
1080 
1081  call h5dread_f(dset_id, h5t_native_integer,dataset,totalshape, hdferr,&
1082  file_space_id = filespace_id, xfer_prp = plist_id, mem_space_id = memspace_id)
1083  if (hdferr < 0) call io_error(1,ext_msg='HDF5_read_int7: h5dread_f')
1084 
1085  call finalize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id)
1086 
1087 end subroutine hdf5_read_int7
1088 
1089 
1090 !--------------------------------------------------------------------------------------------------
1092 !--------------------------------------------------------------------------------------------------
1093 subroutine hdf5_write_real1(loc_id,dataset,datasetName,parallel)
1095  real(pReal), intent(inout), dimension(:) :: dataset
1096  integer(HID_T), intent(in) :: loc_id
1097  character(len=*), intent(in) :: datasetName
1098  logical, intent(in), optional :: parallel
1099 
1100 
1101  integer :: hdferr
1102  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id
1103  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
1104  myStart, &
1105  myShape, & !< shape of the dataset (this process)
1106  totalShape
1107 
1108 !---------------------------------------------------------------------------------------------------
1109 ! determine shape of dataset
1110  myshape = int(shape(dataset),hsize_t)
1111  if (any(myshape(1:size(myshape)-1) == 0)) return
1112 
1113  if (present(parallel)) then
1114  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1115  mystart, totalshape,loc_id,myshape,datasetname,h5t_native_double,parallel)
1116  else
1117  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1118  mystart, totalshape,loc_id,myshape,datasetname,h5t_native_double,.false.)
1119  endif
1120 
1121  if (product(totalshape) /= 0) then
1122  call h5dwrite_f(dset_id, h5t_native_double,dataset,int(totalshape,hsize_t), hdferr,&
1123  file_space_id = filespace_id, mem_space_id = memspace_id, xfer_prp = plist_id)
1124  if (hdferr < 0) call io_error(1,ext_msg='HDF5_write_real1: h5dwrite_f')
1125  endif
1126 
1127  call finalize_write(plist_id, dset_id, filespace_id, memspace_id)
1128 
1129 end subroutine hdf5_write_real1
1130 
1131 !--------------------------------------------------------------------------------------------------
1133 !--------------------------------------------------------------------------------------------------
1134 subroutine hdf5_write_real2(loc_id,dataset,datasetName,parallel)
1136  real(pReal), intent(inout), dimension(:,:) :: dataset
1137  integer(HID_T), intent(in) :: loc_id
1138  character(len=*), intent(in) :: datasetName
1139  logical, intent(in), optional :: parallel
1140 
1141 
1142  integer :: hdferr
1143  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id
1144  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
1145  myStart, &
1146  myShape, & !< shape of the dataset (this process)
1147  totalShape
1148 
1149 !---------------------------------------------------------------------------------------------------
1150 ! determine shape of dataset
1151  myshape = int(shape(dataset),hsize_t)
1152  if (any(myshape(1:size(myshape)-1) == 0)) return
1153 
1154  if (present(parallel)) then
1155  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1156  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_double,parallel)
1157  else
1158  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1159  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_double,.false.)
1160  endif
1161 
1162  if (product(totalshape) /= 0) then
1163  call h5dwrite_f(dset_id, h5t_native_double,dataset,int(totalshape,hsize_t), hdferr,&
1164  file_space_id = filespace_id, mem_space_id = memspace_id, xfer_prp = plist_id)
1165  if (hdferr < 0) call io_error(1,ext_msg='HDF5_write_real2: h5dwrite_f')
1166  endif
1167 
1168  call finalize_write(plist_id, dset_id, filespace_id, memspace_id)
1169 
1170 end subroutine hdf5_write_real2
1171 
1172 !--------------------------------------------------------------------------------------------------
1174 !--------------------------------------------------------------------------------------------------
1175 subroutine hdf5_write_real3(loc_id,dataset,datasetName,parallel)
1177  real(pReal), intent(inout), dimension(:,:,:) :: dataset
1178  integer(HID_T), intent(in) :: loc_id
1179  character(len=*), intent(in) :: datasetName
1180  logical, intent(in), optional :: parallel
1181 
1182 
1183  integer :: hdferr
1184  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id
1185  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
1186  myStart, &
1187  myShape, & !< shape of the dataset (this process)
1188  totalShape
1189 
1190 !---------------------------------------------------------------------------------------------------
1191 ! determine shape of dataset
1192  myshape = int(shape(dataset),hsize_t)
1193  if (any(myshape(1:size(myshape)-1) == 0)) return
1194 
1195  if (present(parallel)) then
1196  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1197  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_double,parallel)
1198  else
1199  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1200  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_double,.false.)
1201  endif
1202 
1203  if (product(totalshape) /= 0) then
1204  call h5dwrite_f(dset_id, h5t_native_double,dataset,int(totalshape,hsize_t), hdferr,&
1205  file_space_id = filespace_id, mem_space_id = memspace_id, xfer_prp = plist_id)
1206  if (hdferr < 0) call io_error(1,ext_msg='HDF5_write_real3: h5dwrite_f')
1207  endif
1208 
1209  call finalize_write(plist_id, dset_id, filespace_id, memspace_id)
1210 
1211 end subroutine hdf5_write_real3
1212 
1213 !--------------------------------------------------------------------------------------------------
1215 !--------------------------------------------------------------------------------------------------
1216 subroutine hdf5_write_real4(loc_id,dataset,datasetName,parallel)
1218  real(pReal), intent(inout), dimension(:,:,:,:) :: dataset
1219  integer(HID_T), intent(in) :: loc_id
1220  character(len=*), intent(in) :: datasetName
1221  logical, intent(in), optional :: parallel
1222 
1223 
1224  integer :: hdferr
1225  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id
1226  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
1227  myStart, &
1228  myShape, & !< shape of the dataset (this process)
1229  totalShape
1230 
1231 !---------------------------------------------------------------------------------------------------
1232 ! determine shape of dataset
1233  myshape = int(shape(dataset),hsize_t)
1234  if (any(myshape(1:size(myshape)-1) == 0)) return
1235 
1236  if (present(parallel)) then
1237  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1238  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_double,parallel)
1239  else
1240  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1241  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_double,.false.)
1242  endif
1243 
1244  if (product(totalshape) /= 0) then
1245  call h5dwrite_f(dset_id, h5t_native_double,dataset,int(totalshape,hsize_t), hdferr,&
1246  file_space_id = filespace_id, mem_space_id = memspace_id, xfer_prp = plist_id)
1247  if (hdferr < 0) call io_error(1,ext_msg='HDF5_write_real4: h5dwrite_f')
1248  endif
1249 
1250  call finalize_write(plist_id, dset_id, filespace_id, memspace_id)
1251 
1252 end subroutine hdf5_write_real4
1253 
1254 
1255 !--------------------------------------------------------------------------------------------------
1257 !--------------------------------------------------------------------------------------------------
1258 subroutine hdf5_write_real5(loc_id,dataset,datasetName,parallel)
1260  real(pReal), intent(inout), dimension(:,:,:,:,:) :: dataset
1261  integer(HID_T), intent(in) :: loc_id
1262  character(len=*), intent(in) :: datasetName
1263  logical, intent(in), optional :: parallel
1264 
1265 
1266  integer :: hdferr
1267  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id
1268  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
1269  myStart, &
1270  myShape, & !< shape of the dataset (this process)
1271  totalShape
1272 
1273 !---------------------------------------------------------------------------------------------------
1274 ! determine shape of dataset
1275  myshape = int(shape(dataset),hsize_t)
1276  if (any(myshape(1:size(myshape)-1) == 0)) return
1277 
1278  if (present(parallel)) then
1279  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1280  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_double,parallel)
1281  else
1282  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1283  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_double,.false.)
1284  endif
1285 
1286  if (product(totalshape) /= 0) then
1287  call h5dwrite_f(dset_id, h5t_native_double,dataset,int(totalshape,hsize_t), hdferr,&
1288  file_space_id = filespace_id, mem_space_id = memspace_id, xfer_prp = plist_id)
1289  if (hdferr < 0) call io_error(1,ext_msg='HDF5_write_real5: h5dwrite_f')
1290  endif
1291 
1292  call finalize_write(plist_id, dset_id, filespace_id, memspace_id)
1293 
1294 end subroutine hdf5_write_real5
1295 
1296 !--------------------------------------------------------------------------------------------------
1298 !--------------------------------------------------------------------------------------------------
1299 subroutine hdf5_write_real6(loc_id,dataset,datasetName,parallel)
1301  real(pReal), intent(inout), dimension(:,:,:,:,:,:) :: dataset
1302  integer(HID_T), intent(in) :: loc_id
1303  character(len=*), intent(in) :: datasetName
1304  logical, intent(in), optional :: parallel
1305 
1306 
1307  integer :: hdferr
1308  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id
1309  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
1310  myStart, &
1311  myShape, & !< shape of the dataset (this process)
1312  totalShape
1313 
1314 !---------------------------------------------------------------------------------------------------
1315 ! determine shape of dataset
1316  myshape = int(shape(dataset),hsize_t)
1317  if (any(myshape(1:size(myshape)-1) == 0)) return
1318 
1319  if (present(parallel)) then
1320  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1321  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_double,parallel)
1322  else
1323  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1324  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_double,.false.)
1325  endif
1326 
1327  if (product(totalshape) /= 0) then
1328  call h5dwrite_f(dset_id, h5t_native_double,dataset,int(totalshape,hsize_t), hdferr,&
1329  file_space_id = filespace_id, mem_space_id = memspace_id, xfer_prp = plist_id)
1330  if (hdferr < 0) call io_error(1,ext_msg='HDF5_write_real6: h5dwrite_f')
1331  endif
1332 
1333  call finalize_write(plist_id, dset_id, filespace_id, memspace_id)
1334 
1335 end subroutine hdf5_write_real6
1336 
1337 !--------------------------------------------------------------------------------------------------
1339 !--------------------------------------------------------------------------------------------------
1340 subroutine hdf5_write_real7(loc_id,dataset,datasetName,parallel)
1342  real(pReal), intent(inout), dimension(:,:,:,:,:,:,:) :: dataset
1343  integer(HID_T), intent(in) :: loc_id
1344  character(len=*), intent(in) :: datasetName
1345  logical, intent(in), optional :: parallel
1346 
1347 
1348  integer :: hdferr
1349  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id
1350  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
1351  myStart, &
1352  myShape, & !< shape of the dataset (this process)
1353  totalShape
1354 
1355 !---------------------------------------------------------------------------------------------------
1356 ! determine shape of dataset
1357  myshape = int(shape(dataset),hsize_t)
1358  if (any(myshape(1:size(myshape)-1) == 0)) return
1359 
1360  if (present(parallel)) then
1361  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1362  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_double,parallel)
1363  else
1364  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1365  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_double,.false.)
1366  endif
1367 
1368  if (product(totalshape) /= 0) then
1369  call h5dwrite_f(dset_id, h5t_native_double,dataset,int(totalshape,hsize_t), hdferr,&
1370  file_space_id = filespace_id, mem_space_id = memspace_id, xfer_prp = plist_id)
1371  if (hdferr < 0) call io_error(1,ext_msg='HDF5_write_real7: h5dwrite_f')
1372  endif
1373 
1374  call finalize_write(plist_id, dset_id, filespace_id, memspace_id)
1375 
1376 end subroutine hdf5_write_real7
1377 
1378 
1379 !--------------------------------------------------------------------------------------------------
1381 !--------------------------------------------------------------------------------------------------
1382 subroutine hdf5_write_int1(loc_id,dataset,datasetName,parallel)
1384  integer, intent(inout), dimension(:) :: dataset
1385  integer(HID_T), intent(in) :: loc_id
1386  character(len=*), intent(in) :: datasetName
1387  logical, intent(in), optional :: parallel
1388 
1389 
1390  integer :: hdferr
1391  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id
1392  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
1393  myStart, &
1394  myShape, & !< shape of the dataset (this process)
1395  totalShape
1396 
1397 !---------------------------------------------------------------------------------------------------
1398 ! determine shape of dataset
1399  myshape = int(shape(dataset),hsize_t)
1400  if (any(myshape(1:size(myshape)-1) == 0)) return
1401 
1402  if (present(parallel)) then
1403  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1404  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_integer,parallel)
1405  else
1406  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1407  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_integer,.false.)
1408  endif
1409 
1410  if (product(totalshape) /= 0) then
1411  call h5dwrite_f(dset_id, h5t_native_integer,dataset,int(totalshape,hsize_t), hdferr,&
1412  file_space_id = filespace_id, mem_space_id = memspace_id, xfer_prp = plist_id)
1413  if (hdferr < 0) call io_error(1,ext_msg='HDF5_write_int1: h5dwrite_f')
1414  endif
1415 
1416  call finalize_write(plist_id, dset_id, filespace_id, memspace_id)
1417 
1418 end subroutine hdf5_write_int1
1419 
1420 !--------------------------------------------------------------------------------------------------
1422 !--------------------------------------------------------------------------------------------------
1423 subroutine hdf5_write_int2(loc_id,dataset,datasetName,parallel)
1425  integer, intent(inout), dimension(:,:) :: dataset
1426  integer(HID_T), intent(in) :: loc_id
1427  character(len=*), intent(in) :: datasetName
1428  logical, intent(in), optional :: parallel
1429 
1430 
1431  integer :: hdferr
1432  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id
1433  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
1434  myStart, &
1435  myShape, & !< shape of the dataset (this process)
1436  totalShape
1437 
1438 !---------------------------------------------------------------------------------------------------
1439 ! determine shape of dataset
1440  myshape = int(shape(dataset),hsize_t)
1441  if (any(myshape(1:size(myshape)-1) == 0)) return
1442 
1443  if (present(parallel)) then
1444  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1445  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_integer,parallel)
1446  else
1447  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1448  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_integer,.false.)
1449  endif
1450 
1451  if (product(totalshape) /= 0) then
1452  call h5dwrite_f(dset_id, h5t_native_integer,dataset,int(totalshape,hsize_t), hdferr,&
1453  file_space_id = filespace_id, mem_space_id = memspace_id, xfer_prp = plist_id)
1454  if (hdferr < 0) call io_error(1,ext_msg='HDF5_write_int2: h5dwrite_f')
1455  endif
1456 
1457  call finalize_write(plist_id, dset_id, filespace_id, memspace_id)
1458 
1459 end subroutine hdf5_write_int2
1460 
1461 !--------------------------------------------------------------------------------------------------
1463 !--------------------------------------------------------------------------------------------------
1464 subroutine hdf5_write_int3(loc_id,dataset,datasetName,parallel)
1466  integer, intent(inout), dimension(:,:,:) :: dataset
1467  integer(HID_T), intent(in) :: loc_id
1468  character(len=*), intent(in) :: datasetName
1469  logical, intent(in), optional :: parallel
1470 
1471 
1472  integer :: hdferr
1473  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id
1474  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
1475  myStart, &
1476  myShape, & !< shape of the dataset (this process)
1477  totalShape
1478 
1479 !---------------------------------------------------------------------------------------------------
1480 ! determine shape of dataset
1481  myshape = int(shape(dataset),hsize_t)
1482  if (any(myshape(1:size(myshape)-1) == 0)) return
1483 
1484  if (present(parallel)) then
1485  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1486  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_integer,parallel)
1487  else
1488  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1489  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_integer,.false.)
1490  endif
1491 
1492  if (product(totalshape) /= 0) then
1493  call h5dwrite_f(dset_id, h5t_native_integer,dataset,int(totalshape,hsize_t), hdferr,&
1494  file_space_id = filespace_id, mem_space_id = memspace_id, xfer_prp = plist_id)
1495  if (hdferr < 0) call io_error(1,ext_msg='HDF5_write_int3: h5dwrite_f')
1496  endif
1497 
1498  call finalize_write(plist_id, dset_id, filespace_id, memspace_id)
1499 
1500 end subroutine hdf5_write_int3
1501 
1502 !--------------------------------------------------------------------------------------------------
1504 !--------------------------------------------------------------------------------------------------
1505 subroutine hdf5_write_int4(loc_id,dataset,datasetName,parallel)
1507  integer, intent(inout), dimension(:,:,:,:) :: dataset
1508  integer(HID_T), intent(in) :: loc_id
1509  character(len=*), intent(in) :: datasetName
1510  logical, intent(in), optional :: parallel
1511 
1512 
1513  integer :: hdferr
1514  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id
1515  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
1516  myStart, &
1517  myShape, & !< shape of the dataset (this process)
1518  totalShape
1519 
1520 !---------------------------------------------------------------------------------------------------
1521 ! determine shape of dataset
1522  myshape = int(shape(dataset),hsize_t)
1523  if (any(myshape(1:size(myshape)-1) == 0)) return
1524 
1525  if (present(parallel)) then
1526  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1527  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_integer,parallel)
1528  else
1529  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1530  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_integer,.false.)
1531  endif
1532 
1533  if (product(totalshape) /= 0) then
1534  call h5dwrite_f(dset_id, h5t_native_integer,dataset,int(totalshape,hsize_t), hdferr,&
1535  file_space_id = filespace_id, mem_space_id = memspace_id, xfer_prp = plist_id)
1536  if (hdferr < 0) call io_error(1,ext_msg='HDF5_write_int4: h5dwrite_f')
1537  endif
1538 
1539  call finalize_write(plist_id, dset_id, filespace_id, memspace_id)
1540 
1541 end subroutine hdf5_write_int4
1542 
1543 !--------------------------------------------------------------------------------------------------
1545 !--------------------------------------------------------------------------------------------------
1546 subroutine hdf5_write_int5(loc_id,dataset,datasetName,parallel)
1548  integer, intent(inout), dimension(:,:,:,:,:) :: dataset
1549  integer(HID_T), intent(in) :: loc_id
1550  character(len=*), intent(in) :: datasetName
1551  logical, intent(in), optional :: parallel
1552 
1553 
1554  integer :: hdferr
1555  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id
1556  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
1557  myStart, &
1558  myShape, & !< shape of the dataset (this process)
1559  totalShape
1560 
1561 !---------------------------------------------------------------------------------------------------
1562 ! determine shape of dataset
1563  myshape = int(shape(dataset),hsize_t)
1564  if (any(myshape(1:size(myshape)-1) == 0)) return
1565 
1566  if (present(parallel)) then
1567  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1568  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_integer,parallel)
1569  else
1570  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1571  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_integer,.false.)
1572  endif
1573 
1574  if (product(totalshape) /= 0) then
1575  call h5dwrite_f(dset_id, h5t_native_integer,dataset,int(totalshape,hsize_t), hdferr,&
1576  file_space_id = filespace_id, mem_space_id = memspace_id, xfer_prp = plist_id)
1577  if (hdferr < 0) call io_error(1,ext_msg='HDF5_write_int5: h5dwrite_f')
1578  endif
1579 
1580  call finalize_write(plist_id, dset_id, filespace_id, memspace_id)
1581 
1582 end subroutine hdf5_write_int5
1583 
1584 !--------------------------------------------------------------------------------------------------
1586 !--------------------------------------------------------------------------------------------------
1587 subroutine hdf5_write_int6(loc_id,dataset,datasetName,parallel)
1589  integer, intent(inout), dimension(:,:,:,:,:,:) :: dataset
1590  integer(HID_T), intent(in) :: loc_id
1591  character(len=*), intent(in) :: datasetName
1592  logical, intent(in), optional :: parallel
1593 
1594 
1595  integer :: hdferr
1596  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id
1597  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
1598  myStart, &
1599  myShape, & !< shape of the dataset (this process)
1600  totalShape
1601 
1602 !---------------------------------------------------------------------------------------------------
1603 ! determine shape of dataset
1604  myshape = int(shape(dataset),hsize_t)
1605  if (any(myshape(1:size(myshape)-1) == 0)) return
1606 
1607  if (present(parallel)) then
1608  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1609  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_integer,parallel)
1610  else
1611  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1612  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_integer,.false.)
1613  endif
1614 
1615  if (product(totalshape) /= 0) then
1616  call h5dwrite_f(dset_id, h5t_native_integer,dataset,int(totalshape,hsize_t), hdferr,&
1617  file_space_id = filespace_id, mem_space_id = memspace_id, xfer_prp = plist_id)
1618  if (hdferr < 0) call io_error(1,ext_msg='HDF5_write_int6: h5dwrite_f')
1619  endif
1620 
1621  call finalize_write(plist_id, dset_id, filespace_id, memspace_id)
1622 
1623 end subroutine hdf5_write_int6
1624 
1625 !--------------------------------------------------------------------------------------------------
1627 !--------------------------------------------------------------------------------------------------
1628 subroutine hdf5_write_int7(loc_id,dataset,datasetName,parallel)
1630  integer, intent(inout), dimension(:,:,:,:,:,:,:) :: dataset
1631  integer(HID_T), intent(in) :: loc_id
1632  character(len=*), intent(in) :: datasetName
1633  logical, intent(in), optional :: parallel
1634 
1635 
1636  integer :: hdferr
1637  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id
1638  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
1639  myStart, &
1640  myShape, & !< shape of the dataset (this process)
1641  totalShape
1642 
1643 !---------------------------------------------------------------------------------------------------
1644 ! determine shape of dataset
1645  myshape = int(shape(dataset),hsize_t)
1646  if (any(myshape(1:size(myshape)-1) == 0)) return
1647 
1648  if (present(parallel)) then
1649  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1650  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_integer,parallel)
1651  else
1652  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1653  mystart, totalshape, loc_id,myshape,datasetname,h5t_native_integer,.false.)
1654  endif
1655 
1656  if (product(totalshape) /= 0) then
1657  call h5dwrite_f(dset_id, h5t_native_integer,dataset,int(totalshape,hsize_t), hdferr,&
1658  file_space_id = filespace_id, mem_space_id = memspace_id, xfer_prp = plist_id)
1659  if (hdferr < 0) call io_error(1,ext_msg='HDF5_write_int7: h5dwrite_f')
1660  endif
1661 
1662  call finalize_write(plist_id, dset_id, filespace_id, memspace_id)
1663 
1664 end subroutine hdf5_write_int7
1665 
1666 
1667 !--------------------------------------------------------------------------------------------------
1669 ! ToDo: It might be possible to write the dataset as a whole
1670 ! ToDo: We could optionally write out other representations (axis angle, euler, ...)
1671 !--------------------------------------------------------------------------------------------------
1672 subroutine hdf5_write_rotation(loc_id,dataset,datasetName,parallel)
1674  type(rotation), intent(in), dimension(:) :: dataset
1675  integer(HID_T), intent(in) :: loc_id
1676  character(len=*), intent(in) :: datasetName
1677  logical, intent(in), optional :: parallel
1678 
1679  integer :: hdferr
1680  real(pReal), dimension(4,size(dataset)) :: dataset_asArray
1681  integer(HID_T) :: dset_id, filespace_id, memspace_id, plist_id,dtype_id,w_id,x_id,y_id,z_id
1682  integer(HSIZE_T), dimension(size(shape(dataset))) :: &
1683  myStart, &
1684  myShape, & !< shape of the dataset (this process)
1685  totalShape
1686  integer(SIZE_T) :: type_size_real
1687  integer :: i
1688 
1689  do i = 1, size(dataset)
1690  dataset_asarray(1:4,i) = dataset(i)%asQuaternion()
1691  enddo
1692 
1693 !---------------------------------------------------------------------------------------------------
1694 ! determine shape of dataset
1695  myshape = int(shape(dataset),hsize_t)
1696 
1697 !---------------------------------------------------------------------------------------------------
1698 ! compound type: name of each quaternion component
1699  call h5tget_size_f(h5t_native_double, type_size_real, hdferr)
1700 
1701  call h5tcreate_f(h5t_compound_f, type_size_real*4_size_t, dtype_id, hdferr)
1702  call h5tinsert_f(dtype_id, "w", type_size_real*0_size_t, h5t_native_double, hdferr)
1703  call h5tinsert_f(dtype_id, "x", type_size_real*1_size_t, h5t_native_double, hdferr)
1704  call h5tinsert_f(dtype_id, "y", type_size_real*2_size_t, h5t_native_double, hdferr)
1705  call h5tinsert_f(dtype_id, "z", type_size_real*3_size_t, h5t_native_double, hdferr)
1706 
1707  if (present(parallel)) then
1708  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1709  mystart, totalshape, loc_id,myshape,datasetname,dtype_id,parallel)
1710  else
1711  call initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1712  mystart, totalshape, loc_id,myshape,datasetname,dtype_id,.false.)
1713  endif
1714 
1715  call h5pset_preserve_f(plist_id, .true., hdferr)
1716 
1717  if (product(totalshape) /= 0) then
1718  call h5tcreate_f(h5t_compound_f, type_size_real, x_id, hdferr)
1719  call h5tinsert_f(x_id, "x", 0_size_t, h5t_native_double, hdferr)
1720  call h5tcreate_f(h5t_compound_f, type_size_real, w_id, hdferr)
1721  call h5tinsert_f(w_id, "w", 0_size_t, h5t_native_double, hdferr)
1722  call h5tcreate_f(h5t_compound_f, type_size_real, y_id, hdferr)
1723  call h5tinsert_f(y_id, "y", 0_size_t, h5t_native_double, hdferr)
1724  call h5tcreate_f(h5t_compound_f, type_size_real, z_id, hdferr)
1725  call h5tinsert_f(z_id, "z", 0_size_t, h5t_native_double, hdferr)
1726 
1727  call h5dwrite_f(dset_id, w_id,dataset_asarray(1,:),int(totalshape,hsize_t), hdferr,&
1728  file_space_id = filespace_id, mem_space_id = memspace_id, xfer_prp = plist_id)
1729  call h5dwrite_f(dset_id, x_id,dataset_asarray(2,:),int(totalshape,hsize_t), hdferr,&
1730  file_space_id = filespace_id, mem_space_id = memspace_id, xfer_prp = plist_id)
1731  call h5dwrite_f(dset_id, y_id,dataset_asarray(3,:),int(totalshape,hsize_t), hdferr,&
1732  file_space_id = filespace_id, mem_space_id = memspace_id, xfer_prp = plist_id)
1733  call h5dwrite_f(dset_id, z_id,dataset_asarray(4,:),int(totalshape,hsize_t), hdferr,&
1734  file_space_id = filespace_id, mem_space_id = memspace_id, xfer_prp = plist_id)
1735  if (hdferr < 0) call io_error(1,ext_msg='HDF5_write_rotation: h5dwrite_f')
1736  endif
1737 
1738  call finalize_write(plist_id, dset_id, filespace_id, memspace_id)
1739 
1740 end subroutine hdf5_write_rotation
1741 
1742 
1743 !--------------------------------------------------------------------------------------------------
1745 !--------------------------------------------------------------------------------------------------
1746 subroutine initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, &
1747  myStart, globalShape, &
1748  loc_id,localShape,datasetName,parallel)
1750  integer(HID_T), intent(in) :: loc_id
1751  character(len=*), intent(in) :: datasetName
1752  logical, intent(in) :: parallel
1753  integer(HSIZE_T), intent(in), dimension(:) :: &
1754  localShape
1755  integer(HSIZE_T), intent(out), dimension(size(localShape,1)):: &
1756  myStart, &
1757  globalShape
1758  integer(HID_T), intent(out) :: dset_id, filespace_id, memspace_id, plist_id, aplist_id
1759 
1760  integer, dimension(worldsize) :: &
1761  readSize
1762  integer :: ierr
1763  integer :: hdferr
1764 
1765 !-------------------------------------------------------------------------------------------------
1766 ! creating a property list for transfer properties (is collective for MPI)
1767  call h5pcreate_f(h5p_dataset_xfer_f, plist_id, hdferr)
1768  if (hdferr < 0) call io_error(1,ext_msg='initialize_read: h5pcreate_f')
1769 
1770 !--------------------------------------------------------------------------------------------------
1771  readsize = 0
1772  readsize(worldrank+1) = int(localshape(ubound(localshape,1)))
1773 
1774  if (parallel) then
1775  call h5pset_dxpl_mpio_f(plist_id, h5fd_mpio_collective_f, hdferr)
1776  if (hdferr < 0) call io_error(1,ext_msg='initialize_read: h5pset_dxpl_mpio_f')
1777  call mpi_allreduce(mpi_in_place,readsize,worldsize,mpi_int,mpi_sum,petsc_comm_world,ierr) ! get total output size over each process
1778  if (ierr /= 0) call io_error(894,ext_msg='initialize_read: MPI_allreduce')
1779  endif
1780 
1781  mystart = int(0,hsize_t)
1782  mystart(ubound(mystart)) = int(sum(readsize(1:worldrank)),hsize_t)
1783  globalshape = [localshape(1:ubound(localshape,1)-1),int(sum(readsize),hsize_t)]
1784 
1785 !--------------------------------------------------------------------------------------------------
1786 ! create dataspace in memory (local shape)
1787  call h5screate_simple_f(size(localshape), localshape, memspace_id, hdferr, localshape)
1788  if (hdferr < 0) call io_error(1,ext_msg='initialize_read: h5screate_simple_f/memspace_id')
1789 
1790 !--------------------------------------------------------------------------------------------------
1791 ! creating a property list for IO and set it to collective
1792  call h5pcreate_f(h5p_dataset_access_f, aplist_id, hdferr)
1793  if (hdferr < 0) call io_error(1,ext_msg='initialize_read: h5pcreate_f')
1794 
1795  call h5pset_all_coll_metadata_ops_f(aplist_id, .true., hdferr)
1796  if (hdferr < 0) call io_error(1,ext_msg='initialize_read: h5pset_all_coll_metadata_ops_f')
1797 
1798 
1799 !--------------------------------------------------------------------------------------------------
1800 ! open the dataset in the file and get the space ID
1801  call h5dopen_f(loc_id,datasetname,dset_id,hdferr, dapl_id = aplist_id)
1802  if (hdferr < 0) call io_error(1,ext_msg='initialize_read: h5dopen_f')
1803  call h5dget_space_f(dset_id, filespace_id, hdferr)
1804  if (hdferr < 0) call io_error(1,ext_msg='initialize_read: h5dget_space_f')
1805 
1806 !--------------------------------------------------------------------------------------------------
1807 ! select a hyperslab (the portion of the current process) in the file
1808  call h5sselect_hyperslab_f(filespace_id, h5s_select_set_f, mystart, localshape, hdferr)
1809  if (hdferr < 0) call io_error(1,ext_msg='initialize_read: h5sselect_hyperslab_f')
1810 
1811 end subroutine initialize_read
1812 
1813 
1814 !--------------------------------------------------------------------------------------------------
1816 !--------------------------------------------------------------------------------------------------
1817 subroutine finalize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id)
1819  integer(HID_T), intent(in) :: dset_id, filespace_id, memspace_id, plist_id, aplist_id
1820  integer :: hdferr
1821 
1822  call h5pclose_f(plist_id, hdferr)
1823  if (hdferr < 0) call io_error(1,ext_msg='finalize_read: plist_id')
1824  call h5pclose_f(aplist_id, hdferr)
1825  if (hdferr < 0) call io_error(1,ext_msg='finalize_read: aplist_id')
1826  call h5dclose_f(dset_id, hdferr)
1827  if (hdferr < 0) call io_error(1,ext_msg='finalize_read: h5dclose_f')
1828  call h5sclose_f(filespace_id, hdferr)
1829  if (hdferr < 0) call io_error(1,ext_msg='finalize_read: h5sclose_f/filespace_id')
1830  call h5sclose_f(memspace_id, hdferr)
1831  if (hdferr < 0) call io_error(1,ext_msg='finalize_read: h5sclose_f/memspace_id')
1832 
1833 end subroutine finalize_read
1834 
1835 
1836 !--------------------------------------------------------------------------------------------------
1838 !--------------------------------------------------------------------------------------------------
1839 subroutine initialize_write(dset_id, filespace_id, memspace_id, plist_id, &
1840  myStart, totalShape, &
1841  loc_id,myShape,datasetName,datatype,parallel)
1843  integer(HID_T), intent(in) :: loc_id
1844  character(len=*), intent(in) :: datasetName
1845  logical, intent(in) :: parallel
1846  integer(HID_T), intent(in) :: datatype
1847  integer(HSIZE_T), intent(in), dimension(:) :: &
1848  myShape
1849  integer(HSIZE_T), intent(out), dimension(size(myShape,1)):: &
1850  myStart, &
1851  totalShape
1852  integer(HID_T), intent(out) :: dset_id, filespace_id, memspace_id, plist_id
1853 
1854  integer, dimension(worldsize) :: &
1855  writeSize
1856  integer :: ierr
1857  integer :: hdferr
1858 
1859 !-------------------------------------------------------------------------------------------------
1860 ! creating a property list for transfer properties (is collective when reading in parallel)
1861  call h5pcreate_f(h5p_dataset_xfer_f, plist_id, hdferr)
1862  if (hdferr < 0) call io_error(1,ext_msg='initialize_write: h5pcreate_f')
1863 
1864  if (parallel) then
1865  call h5pset_dxpl_mpio_f(plist_id, h5fd_mpio_collective_f, hdferr)
1866  if (hdferr < 0) call io_error(1,ext_msg='initialize_write: h5pset_dxpl_mpio_f')
1867  endif
1868 
1869 
1870 !--------------------------------------------------------------------------------------------------
1871 ! determine the global data layout among all processes
1872  writesize = 0
1873  writesize(worldrank+1) = int(myshape(ubound(myshape,1)))
1874 
1875  if (parallel) then
1876  call mpi_allreduce(mpi_in_place,writesize,worldsize,mpi_int,mpi_sum,petsc_comm_world,ierr) ! get total output size over each process
1877  if (ierr /= 0) call io_error(894,ext_msg='initialize_write: MPI_allreduce')
1878  endif
1879 
1880  mystart = int(0,hsize_t)
1881  mystart(ubound(mystart)) = int(sum(writesize(1:worldrank)),hsize_t)
1882  totalshape = [myshape(1:ubound(myshape,1)-1),int(sum(writesize),hsize_t)]
1883 
1884 !--------------------------------------------------------------------------------------------------
1885 ! create dataspace in memory (local shape) and in file (global shape)
1886  call h5screate_simple_f(size(myshape), myshape, memspace_id, hdferr, myshape)
1887  if (hdferr < 0) call io_error(1,ext_msg='initialize_write: h5dopen_f')
1888  call h5screate_simple_f(size(totalshape), totalshape, filespace_id, hdferr, totalshape)
1889  if (hdferr < 0) call io_error(1,ext_msg='initialize_write: h5dget_space_f')
1890 
1891 !--------------------------------------------------------------------------------------------------
1892 ! create dataset in the file and select a hyperslab from it (the portion of the current process)
1893  call h5dcreate_f(loc_id, trim(datasetname), datatype, filespace_id, dset_id, hdferr)
1894  if (hdferr < 0) call io_error(1,ext_msg='initialize_write: h5dcreate_f')
1895  call h5sselect_hyperslab_f(filespace_id, h5s_select_set_f, mystart, myshape, hdferr)
1896  if (hdferr < 0) call io_error(1,ext_msg='initialize_write: h5sselect_hyperslab_f')
1897 
1898 end subroutine initialize_write
1899 
1900 
1901 !--------------------------------------------------------------------------------------------------
1903 !--------------------------------------------------------------------------------------------------
1904 subroutine finalize_write(plist_id, dset_id, filespace_id, memspace_id)
1906  integer(HID_T), intent(in) :: dset_id, filespace_id, memspace_id, plist_id
1907  integer :: hdferr
1908 
1909  call h5pclose_f(plist_id, hdferr)
1910  if (hdferr < 0) call io_error(1,ext_msg='finalize_write: plist_id')
1911  call h5dclose_f(dset_id, hdferr)
1912  if (hdferr < 0) call io_error(1,ext_msg='finalize_write: h5dclose_f')
1913  call h5sclose_f(filespace_id, hdferr)
1914  if (hdferr < 0) call io_error(1,ext_msg='finalize_write: h5sclose_f/filespace_id')
1915  call h5sclose_f(memspace_id, hdferr)
1916  if (hdferr < 0) call io_error(1,ext_msg='finalize_write: h5sclose_f/memspace_id')
1917 
1918 end subroutine finalize_write
1919 
1920 end module hdf5_utilities
hdf5_utilities::hdf5_write_real1
subroutine hdf5_write_real1(loc_id, dataset, datasetName, parallel)
write dataset of type real with 1 dimension
Definition: HDF5_utilities.f90:1094
hdf5_utilities::hdf5_addgroup
integer(hid_t) function hdf5_addgroup(fileHandle, groupName)
adds a new group to the fileHandle
Definition: HDF5_utilities.f90:182
hdf5_utilities::hdf5_read_int4
subroutine hdf5_read_int4(loc_id, dataset, datasetName, parallel)
read dataset of type integer withh 4 dimensions
Definition: HDF5_utilities.f90:933
hdf5_utilities::hdf5_write_int7
subroutine hdf5_write_int7(loc_id, dataset, datasetName, parallel)
write dataset of type integer with 7 dimensions
Definition: HDF5_utilities.f90:1629
hdf5_utilities::hdf5_read_int7
subroutine hdf5_read_int7(loc_id, dataset, datasetName, parallel)
read dataset of type integer with 7 dimensions
Definition: HDF5_utilities.f90:1053
rotations
rotation storage and conversion
Definition: rotations.f90:53
hdf5_utilities::hdf5_addattribute_str
subroutine hdf5_addattribute_str(loc_id, attrLabel, attrValue, path)
adds a string attribute to the path given relative to the location
Definition: HDF5_utilities.f90:293
io::io_error
subroutine, public io_error(error_ID, el, ip, g, instance, ext_msg)
write error statements to standard out and terminate the Marc/spectral run with exit #9xxx
Definition: IO.f90:305
hdf5_utilities::hdf5_closegroup
subroutine hdf5_closegroup(group_id)
close a group
Definition: HDF5_utilities.f90:251
hdf5_utilities::hdf5_addattribute
attached attributes of type char, integer or real to a file/dataset/group
Definition: HDF5_utilities.f90:76
rotations::rotation
Definition: rotations.f90:63
hdf5_utilities::finalize_read
subroutine finalize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id)
closes HDF5 handles
Definition: HDF5_utilities.f90:1818
hdf5_utilities::hdf5_read_int1
subroutine hdf5_read_int1(loc_id, dataset, datasetName, parallel)
read dataset of type integer with 1 dimension
Definition: HDF5_utilities.f90:812
hdf5_utilities::hdf5_write_real3
subroutine hdf5_write_real3(loc_id, dataset, datasetName, parallel)
write dataset of type real with 3 dimensions
Definition: HDF5_utilities.f90:1176
hdf5_utilities
Definition: HDF5_utilities.f90:11
hdf5_utilities::hdf5_write_int6
subroutine hdf5_write_int6(loc_id, dataset, datasetName, parallel)
write dataset of type integer with 6 dimensions
Definition: HDF5_utilities.f90:1588
hdf5_utilities::hdf5_openfile
integer(hid_t) function hdf5_openfile(fileName, mode, parallel)
open and initializes HDF5 output file
Definition: HDF5_utilities.f90:119
numerics::worldsize
integer, public, protected worldsize
MPI worldsize (/=1 for MPI simulations only)
Definition: numerics.f90:1470
hdf5_utilities::hdf5_read_int2
subroutine hdf5_read_int2(loc_id, dataset, datasetName, parallel)
read dataset of type integer with 2 dimensions
Definition: HDF5_utilities.f90:853
hdf5_utilities::initialize_write
subroutine initialize_write(dset_id, filespace_id, memspace_id, plist_id, myStart, totalShape, loc_id, myShape, datasetName, datatype, parallel)
initialize HDF5 handles, determines global shape and start for parallel write
Definition: HDF5_utilities.f90:1842
hdf5_utilities::hdf5_addattribute_int_array
subroutine hdf5_addattribute_int_array(loc_id, attrLabel, attrValue, path)
adds a integer attribute to the path given relative to the location
Definition: HDF5_utilities.f90:421
hdf5_utilities::initialize_read
subroutine initialize_read(dset_id, filespace_id, memspace_id, plist_id, aplist_id, myStart, globalShape, loc_id, localShape, datasetName, parallel)
initialize HDF5 handles, determines global shape and start for parallel read
Definition: HDF5_utilities.f90:1749
hdf5_utilities::hdf5_write_int5
subroutine hdf5_write_int5(loc_id, dataset, datasetName, parallel)
write dataset of type integer with 5 dimensions
Definition: HDF5_utilities.f90:1547
prec
setting precision for real and int type
Definition: prec.f90:13
hdf5_utilities::hdf5_write_real5
subroutine hdf5_write_real5(loc_id, dataset, datasetName, parallel)
write dataset of type real with 5 dimensions
Definition: HDF5_utilities.f90:1259
hdf5_utilities::hdf5_read_real3
subroutine hdf5_read_real3(loc_id, dataset, datasetName, parallel)
read dataset of type real with 2 dimensions
Definition: HDF5_utilities.f90:611
hdf5_utilities::hdf5_objectexists
logical function hdf5_objectexists(loc_id, path)
check whether a group or a dataset exists
Definition: HDF5_utilities.f90:265
hdf5_utilities::hdf5_read
reads integer or float data of defined shape from file ! ToDo: order of arguments wrong
Definition: HDF5_utilities.f90:29
io
input/output functions, partly depending on chosen solver
Definition: IO.f90:12
hdf5_utilities::hdf5_read_real6
subroutine hdf5_read_real6(loc_id, dataset, datasetName, parallel)
read dataset of type real with 6 dimensions
Definition: HDF5_utilities.f90:731
hdf5_utilities::hdf5_read_real4
subroutine hdf5_read_real4(loc_id, dataset, datasetName, parallel)
read dataset of type real with 4 dimensions
Definition: HDF5_utilities.f90:651
hdf5_utilities::hdf5_read_real2
subroutine hdf5_read_real2(loc_id, dataset, datasetName, parallel)
read dataset of type real with 2 dimensions
Definition: HDF5_utilities.f90:571
hdf5_utilities::hdf5_utilities_init
subroutine hdf5_utilities_init
open libary and do sanity checks
Definition: HDF5_utilities.f90:91
hdf5_utilities::hdf5_write_int3
subroutine hdf5_write_int3(loc_id, dataset, datasetName, parallel)
write dataset of type integer with 3 dimensions
Definition: HDF5_utilities.f90:1465
hdf5_utilities::hdf5_read_int6
subroutine hdf5_read_int6(loc_id, dataset, datasetName, parallel)
read dataset of type integer with 6 dimensions
Definition: HDF5_utilities.f90:1013
hdf5_utilities::hdf5_write_int1
subroutine hdf5_write_int1(loc_id, dataset, datasetName, parallel)
write dataset of type integer with 1 dimension
Definition: HDF5_utilities.f90:1383
hdf5_utilities::hdf5_addattribute_real
subroutine hdf5_addattribute_real(loc_id, attrLabel, attrValue, path)
adds a integer attribute to the path given relative to the location
Definition: HDF5_utilities.f90:380
hdf5_utilities::hdf5_addattribute_real_array
subroutine hdf5_addattribute_real_array(loc_id, attrLabel, attrValue, path)
adds a real attribute to the path given relative to the location
Definition: HDF5_utilities.f90:465
hdf5_utilities::finalize_write
subroutine finalize_write(plist_id, dset_id, filespace_id, memspace_id)
closes HDF5 handles
Definition: HDF5_utilities.f90:1905
hdf5_utilities::hdf5_read_real1
subroutine hdf5_read_real1(loc_id, dataset, datasetName, parallel)
read dataset of type real with 1 dimension
Definition: HDF5_utilities.f90:531
hdf5_utilities::hdf5_read_int3
subroutine hdf5_read_int3(loc_id, dataset, datasetName, parallel)
read dataset of type integer with 3 dimensions
Definition: HDF5_utilities.f90:893
hdf5_utilities::hdf5_write_int2
subroutine hdf5_write_int2(loc_id, dataset, datasetName, parallel)
write dataset of type integer with 2 dimensions
Definition: HDF5_utilities.f90:1424
hdf5_utilities::hdf5_closefile
subroutine hdf5_closefile(fileHandle)
close the opened HDF5 output file
Definition: HDF5_utilities.f90:167
numerics::worldrank
integer, public, protected worldrank
MPI worldrank (/=0 for MPI simulations only)
Definition: numerics.f90:1470
hdf5_utilities::hdf5_write_real6
subroutine hdf5_write_real6(loc_id, dataset, datasetName, parallel)
write dataset of type real with 6 dimensions
Definition: HDF5_utilities.f90:1300
hdf5_utilities::hdf5_write_rotation
subroutine hdf5_write_rotation(loc_id, dataset, datasetName, parallel)
writes a scalar orientation dataset
Definition: HDF5_utilities.f90:1673
hdf5_utilities::hdf5_read_real5
subroutine hdf5_read_real5(loc_id, dataset, datasetName, parallel)
read dataset of type real with 5 dimensions
Definition: HDF5_utilities.f90:691
hdf5_utilities::hdf5_write_real4
subroutine hdf5_write_real4(loc_id, dataset, datasetName, parallel)
write dataset of type real with 4 dimensions
Definition: HDF5_utilities.f90:1217
hdf5_utilities::hdf5_setlink
subroutine hdf5_setlink(loc_id, target_name, link_name)
set link to object in results file
Definition: HDF5_utilities.f90:509
hdf5_utilities::hdf5_write_real2
subroutine hdf5_write_real2(loc_id, dataset, datasetName, parallel)
write dataset of type real with 2 dimensions
Definition: HDF5_utilities.f90:1135
hdf5_utilities::hdf5_read_int5
subroutine hdf5_read_int5(loc_id, dataset, datasetName, parallel)
read dataset of type integer with 5 dimensions
Definition: HDF5_utilities.f90:973
numerics
Managing of parameters related to numerics.
Definition: numerics.f90:10
hdf5_utilities::hdf5_opengroup
integer(hid_t) function hdf5_opengroup(fileHandle, groupName)
open an existing group of a file
Definition: HDF5_utilities.f90:215
hdf5_utilities::hdf5_write_real7
subroutine hdf5_write_real7(loc_id, dataset, datasetName, parallel)
write dataset of type real with 7 dimensions
Definition: HDF5_utilities.f90:1341
hdf5_utilities::hdf5_write_int4
subroutine hdf5_write_int4(loc_id, dataset, datasetName, parallel)
write dataset of type integer with 4 dimensions
Definition: HDF5_utilities.f90:1506
hdf5_utilities::hdf5_addattribute_int
subroutine hdf5_addattribute_int(loc_id, attrLabel, attrValue, path)
adds a integer attribute to the path given relative to the location
Definition: HDF5_utilities.f90:339
hdf5_utilities::hdf5_read_real7
subroutine hdf5_read_real7(loc_id, dataset, datasetName, parallel)
read dataset of type real with 7 dimensions
Definition: HDF5_utilities.f90:771
hdf5_utilities::hdf5_write
writes integer or real data of defined shape to file ! ToDo: order of arguments wrong
Definition: HDF5_utilities.f90:52