Neko  0.8.99
A portable framework for high-order spectral element flow simulations
hdf5_file.F90
Go to the documentation of this file.
1 ! Copyright (c) 2024, The Neko Authors
2 ! All rights reserved.
3 !
4 ! Redistribution and use in source and binary forms, with or without
5 ! modification, are permitted provided that the following conditions
6 ! are met:
7 !
8 ! * Redistributions of source code must retain the above copyright
9 ! notice, this list of conditions and the following disclaimer.
10 !
11 ! * Redistributions in binary form must reproduce the above
12 ! copyright notice, this list of conditions and the following
13 ! disclaimer in the documentation and/or other materials provided
14 ! with the distribution.
15 !
16 ! * Neither the name of the authors nor the names of its
17 ! contributors may be used to endorse or promote products derived
18 ! from this software without specific prior written permission.
19 !
20 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24 ! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25 ! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26 ! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27 ! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28 ! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29 ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30 ! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31 ! POSSIBILITY OF SUCH DAMAGE.
32 !
34 module hdf5_file
35  use num_types, only : rp
36  use generic_file, only : generic_file_t
37  use checkpoint, only : chkp_t
39  use mesh, only : mesh_t
40  use field, only : field_t, field_ptr_t
41  use field_list, only : field_list_t
43  use dofmap, only : dofmap_t
44  use logger
45  use comm
46  use mpi_f08
47 #ifdef HAVE_HDF5
48  use hdf5
49 #endif
50  implicit none
51  private
52 
54  type, public, extends(generic_file_t) :: hdf5_file_t
55  contains
56  procedure :: read => hdf5_file_read
57  procedure :: write => hdf5_file_write
58  end type hdf5_file_t
59 
60 contains
61 
62 #ifdef HAVE_HDF5
63 
65  subroutine hdf5_file_write(this, data, t)
66  class(hdf5_file_t), intent(inout) :: this
67  class(*), target, intent(in) :: data
68  real(kind=rp), intent(in), optional :: t
69  type(mesh_t), pointer :: msh
70  type(dofmap_t), pointer :: dof
71  type(field_ptr_t), allocatable :: fp(:)
72  type(field_series_ptr_t), allocatable :: fsp(:)
73  real(kind=rp), pointer :: dtlag(:)
74  real(kind=rp), pointer :: tlag(:)
75  integer :: ierr, info, drank, i, j
76  integer(hid_t) :: plist_id, file_id, dset_id, grp_id, attr_id
77  integer(hid_t) :: filespace, memspace
78  integer(hsize_t), dimension(1) :: ddim, dcount, doffset
79  integer :: suffix_pos
80  character(len=5) :: id_str
81  character(len=1024) :: fname
82 
83  call hdf5_file_determine_data(data, msh, dof, fp, fsp, dtlag, tlag)
84 
85  suffix_pos = filename_suffix_pos(this%fname)
86  write(id_str, '(i5.5)') this%counter
87  fname = trim(this%fname(1:suffix_pos-1))//id_str//'.h5'
88 
89  call h5open_f(ierr)
90  call h5pcreate_f(h5p_file_access_f, plist_id, ierr)
91  info = mpi_info_null%mpi_val
92  call h5pset_fapl_mpio_f(plist_id, neko_comm%mpi_val, info, ierr)
93 
94  call h5fcreate_f(fname, h5f_acc_trunc_f, &
95  file_id, ierr, access_prp = plist_id)
96 
97  call h5pcreate_f(h5p_dataset_xfer_f, plist_id, ierr)
98  call h5pset_dxpl_mpio_f(plist_id, h5fd_mpio_collective_f, ierr)
99 
100  call h5screate_f(h5s_scalar_f, filespace, ierr)
101  ddim = 1
102 
103  if (present(t)) then
104  call h5acreate_f(file_id, "Time", h5t_native_double, filespace, attr_id, &
105  ierr, h5p_default_f, h5p_default_f)
106  call h5awrite_f(attr_id, h5t_native_double, t, ddim, ierr)
107  call h5aclose_f(attr_id, ierr)
108  end if
109 
110  if (associated(dof)) then
111  call h5acreate_f(file_id, "Lx", h5t_native_integer, filespace, attr_id, &
112  ierr, h5p_default_f, h5p_default_f)
113  call h5awrite_f(attr_id, h5t_native_integer, dof%Xh%lx, ddim, ierr)
114  call h5aclose_f(attr_id, ierr)
115  end if
116 
117  if (associated(msh)) then
118  call h5gcreate_f(file_id, "Mesh", grp_id, ierr, lcpl_id=h5p_default_f, &
119  gcpl_id=h5p_default_f, gapl_id=h5p_default_f)
120 
121  call h5acreate_f(grp_id, "Elements", h5t_native_integer, filespace, attr_id, &
122  ierr, h5p_default_f, h5p_default_f)
123  call h5awrite_f(attr_id, h5t_native_integer, msh%glb_nelv, ddim, ierr)
124  call h5aclose_f(attr_id, ierr)
125 
126  call h5acreate_f(grp_id, "Dimension", h5t_native_integer, filespace, attr_id, &
127  ierr, h5p_default_f, h5p_default_f)
128  call h5awrite_f(attr_id, h5t_native_integer, msh%gdim, ddim, ierr)
129  call h5aclose_f(attr_id, ierr)
130 
131  call h5gclose_f(grp_id, ierr)
132  end if
133 
134 
135  call h5sclose_f(filespace, ierr)
136 
137  !
138  ! Write restart group (tlag, dtlag)
139  !
140  if (associated(tlag) .and. associated(dtlag)) then
141  call h5gcreate_f(file_id, "Restart", grp_id, ierr, lcpl_id=h5p_default_f, &
142  gcpl_id=h5p_default_f, gapl_id=h5p_default_f)
143 
144  drank = 1
145  ddim = size(tlag)
146  doffset(1) = 0
147  if (pe_rank .eq. 0) then
148  dcount = size(tlag)
149  else
150  dcount = 0
151  end if
152 
153  call h5screate_simple_f(drank, ddim, filespace, ierr)
154 
155  call h5dcreate_f(grp_id,'tlag', h5t_native_double, &
156  filespace, dset_id, ierr)
157  call h5dget_space_f(dset_id, filespace, ierr)
158  call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
159  doffset, dcount, ierr)
160  call h5dwrite_f(dset_id, h5t_native_double, tlag, &
161  ddim, ierr, xfer_prp = plist_id)
162  call h5dclose_f(dset_id, ierr)
163 
164  call h5dcreate_f(grp_id,'dtlag', h5t_native_double, &
165  filespace, dset_id, ierr)
166  call h5dget_space_f(dset_id, filespace, ierr)
167  call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
168  doffset, dcount, ierr)
169  call h5dwrite_f(dset_id, h5t_native_double, dtlag, &
170  ddim, ierr, xfer_prp = plist_id)
171  call h5dclose_f(dset_id, ierr)
172 
173  call h5sclose_f(filespace, ierr)
174  call h5gclose_f(grp_id, ierr)
175 
176  end if
177 
178 
179  !
180  ! Write fields group
181  !
182  if (allocated(fp) .or. allocated(fsp)) then
183  call h5gcreate_f(file_id, "Fields", grp_id, ierr, lcpl_id=h5p_default_f, &
184  gcpl_id=h5p_default_f, gapl_id=h5p_default_f)
185 
186  dcount(1) = int(dof%size(), 8)
187  doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
188  ddim = int(dof%size(), 8)
189  drank = 1
190  call mpi_allreduce(mpi_in_place, ddim(1), 1, &
191  mpi_integer8, mpi_sum, neko_comm, ierr)
192 
193  call h5screate_simple_f(drank, ddim, filespace, ierr)
194  call h5screate_simple_f(drank, dcount, memspace, ierr)
195 
196 
197  if (allocated(fp)) then
198  do i = 1, size(fp)
199  call h5dcreate_f(grp_id, fp(i)%ptr%name, h5t_native_double, &
200  filespace, dset_id, ierr)
201  call h5dget_space_f(dset_id, filespace, ierr)
202  call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
203  doffset, dcount, ierr)
204  call h5dwrite_f(dset_id, h5t_native_double, &
205  fp(i)%ptr%x(1,1,1,1), &
206  ddim, ierr, file_space_id = filespace, &
207  mem_space_id = memspace, xfer_prp = plist_id)
208  call h5dclose_f(dset_id, ierr)
209  end do
210  deallocate(fp)
211  end if
212 
213  if (allocated(fsp)) then
214  do i = 1, size(fsp)
215  do j = 1, fsp(i)%ptr%size()
216  call h5dcreate_f(grp_id, fsp(i)%ptr%lf(j)%name, &
217  h5t_native_double, filespace, dset_id, ierr)
218  call h5dget_space_f(dset_id, filespace, ierr)
219  call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
220  doffset, dcount, ierr)
221  call h5dwrite_f(dset_id, h5t_native_double, &
222  fsp(i)%ptr%lf(j)%x(1,1,1,1), &
223  ddim, ierr, file_space_id = filespace, &
224  mem_space_id = memspace, xfer_prp = plist_id)
225  call h5dclose_f(dset_id, ierr)
226  end do
227  end do
228  deallocate(fsp)
229  end if
230 
231  call h5gclose_f(grp_id, ierr)
232  call h5sclose_f(filespace, ierr)
233  call h5sclose_f(memspace, ierr)
234  end if
235 
236  call h5pclose_f(plist_id, ierr)
237  call h5fclose_f(file_id, ierr)
238 
239  call h5close_f(ierr)
240 
241  this%counter = this%counter + 1
242 
243  end subroutine hdf5_file_write
244 
246  subroutine hdf5_file_read(this, data)
247  class(hdf5_file_t) :: this
248  class(*), target, intent(inout) :: data
249  integer(hid_t) :: plist_id, file_id, dset_id, grp_id, attr_id
250  integer(hid_t) :: filespace, memspace
251  integer(hsize_t), dimension(1) :: ddim, dcount, doffset
252  integer :: i,j, ierr, info, glb_nelv, gdim, lx, drank
253  type(mesh_t), pointer :: msh
254  type(dofmap_t), pointer :: dof
255  type(field_ptr_t), allocatable :: fp(:)
256  type(field_series_ptr_t), allocatable :: fsp(:)
257  real(kind=rp), pointer :: dtlag(:)
258  real(kind=rp), pointer :: tlag(:)
259  real(kind=rp) :: t
260 
261  call hdf5_file_determine_data(data, msh, dof, fp, fsp, dtlag, tlag)
262 
263  call h5open_f(ierr)
264  call h5pcreate_f(h5p_file_access_f, plist_id, ierr)
265  info = mpi_info_null%mpi_val
266  call h5pset_fapl_mpio_f(plist_id, neko_comm%mpi_val, info, ierr)
267 
268  call h5fopen_f(trim(this%fname), h5f_acc_rdonly_f, &
269  file_id, ierr, access_prp = plist_id)
270 
271  call h5pcreate_f(h5p_dataset_xfer_f, plist_id, ierr)
272  call h5pset_dxpl_mpio_f(plist_id, h5fd_mpio_collective_f, ierr)
273 
274  ddim = 1
275  call h5aopen_name_f(file_id, 'Time', attr_id, ierr)
276  call h5aread_f(attr_id, h5t_native_double, t, ddim, ierr)
277  call h5aclose_f(attr_id, ierr)
278 
279  select type(data)
280  type is(chkp_t)
281  data%t = t
282  end select
283 
284  call h5aopen_name_f(file_id, 'Lx', attr_id, ierr)
285  call h5aread_f(attr_id, h5t_native_integer, lx, ddim, ierr)
286  call h5aclose_f(attr_id, ierr)
287 
288  call h5gopen_f(file_id, 'Mesh', grp_id, ierr, gapl_id=h5p_default_f)
289 
290  call h5aopen_name_f(grp_id, 'Elements', attr_id, ierr)
291  call h5aread_f(attr_id, h5t_native_integer, glb_nelv, ddim, ierr)
292  call h5aclose_f(attr_id, ierr)
293 
294  call h5aopen_name_f(grp_id, 'Dimension', attr_id, ierr)
295  call h5aread_f(attr_id, h5t_native_integer, gdim, ddim, ierr)
296  call h5aclose_f(attr_id, ierr)
297  call h5gclose_f(grp_id, ierr)
298 
299 
300  if (associated(tlag) .and. associated(dtlag)) then
301  drank = 1
302  ddim = size(tlag)
303  doffset(1) = 0
304  if (pe_rank .eq. 0) then
305  dcount = size(tlag)
306  else
307  dcount = 0
308  end if
309 
310  call h5gopen_f(file_id, 'Restart', grp_id, ierr, gapl_id=h5p_default_f)
311  call h5dopen_f(grp_id, 'tlag', dset_id, ierr)
312  call h5dget_space_f(dset_id, filespace, ierr)
313  call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
314  doffset, dcount, ierr)
315  call h5dread_f(dset_id, h5t_native_double, tlag, ddim, ierr, xfer_prp=plist_id)
316  call h5dclose_f(dset_id, ierr)
317  call h5sclose_f(filespace, ierr)
318 
319  call h5dopen_f(grp_id, 'dtlag', dset_id, ierr)
320  call h5dget_space_f(dset_id, filespace, ierr)
321  call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
322  doffset, dcount, ierr)
323  call h5dread_f(dset_id, h5t_native_double, dtlag, ddim, ierr, xfer_prp=plist_id)
324  call h5dclose_f(dset_id, ierr)
325  call h5sclose_f(filespace, ierr)
326 
327  call h5gclose_f(grp_id, ierr)
328  end if
329 
330  if (allocated(fp) .or. allocated(fsp)) then
331  call h5gopen_f(file_id, 'Fields', grp_id, ierr, gapl_id=h5p_default_f)
332 
333  dcount(1) = int(dof%size(), 8)
334  doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
335  ddim = int(dof%size(), 8)
336  drank = 1
337 
338  dcount(1) = int(dof%size(), 8)
339  doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
340  ddim = int(dof%size(), 8)
341  drank = 1
342  call mpi_allreduce(mpi_in_place, ddim(1), 1, &
343  mpi_integer8, mpi_sum, neko_comm, ierr)
344 
345  call h5screate_simple_f(drank, dcount, memspace, ierr)
346 
347  if (allocated(fp)) then
348  do i = 1, size(fp)
349  call h5dopen_f(grp_id, fp(i)%ptr%name, dset_id, ierr)
350  call h5dget_space_f(dset_id, filespace, ierr)
351  call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
352  doffset, dcount, ierr)
353  call h5dread_f(dset_id, h5t_native_double, &
354  fp(i)%ptr%x(1,1,1,1), &
355  ddim, ierr, file_space_id = filespace, &
356  mem_space_id = memspace, xfer_prp=plist_id)
357  call h5dclose_f(dset_id, ierr)
358  call h5sclose_f(filespace, ierr)
359  end do
360  end if
361 
362  if (allocated(fsp)) then
363  do i = 1, size(fsp)
364  do j = 1, fsp(i)%ptr%size()
365  call h5dopen_f(grp_id, fsp(i)%ptr%lf(j)%name, dset_id, ierr)
366  call h5dget_space_f(dset_id, filespace, ierr)
367  call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
368  doffset, dcount, ierr)
369  call h5dread_f(dset_id, h5t_native_double, &
370  fsp(i)%ptr%lf(j)%x(1,1,1,1), &
371  ddim, ierr, file_space_id = filespace, &
372  mem_space_id = memspace, xfer_prp=plist_id)
373  call h5dclose_f(dset_id, ierr)
374  call h5sclose_f(filespace, ierr)
375  end do
376  end do
377  end if
378  call h5sclose_f(memspace, ierr)
379  call h5gclose_f(grp_id, ierr)
380  end if
381 
382  call h5pclose_f(plist_id, ierr)
383  call h5fclose_f(file_id, ierr)
384 
385  call h5close_f(ierr)
386 
387  end subroutine hdf5_file_read
388 
389 
390  subroutine hdf5_file_determine_data(data, msh, dof, fp, fsp, dtlag, tlag)
391  class(*), target, intent(in) :: data
392  type(mesh_t), pointer, intent(inout) :: msh
393  type(dofmap_t), pointer, intent(inout) :: dof
394  type(field_ptr_t), allocatable, intent(inout) :: fp(:)
395  type(field_series_ptr_t), allocatable, intent(inout) :: fsp(:)
396  real(kind=rp), pointer, intent(inout) :: dtlag(:)
397  real(kind=rp), pointer, intent(inout) :: tlag(:)
398  integer :: i, j, fp_size, fp_cur, fsp_size, fsp_cur
399 
400  select type(data)
401  type is (field_t)
402  dof => data%dof
403  msh => data%msh
404  fp_size = 1
405  allocate(fp(fp_size))
406  fp(1)%ptr => data
407 
408  nullify(dtlag)
409  nullify(tlag)
410 
411  type is (field_list_t)
412 
413  if (data%size() .gt. 0) then
414  allocate(fp(data%size()))
415 
416  dof => data%dof(1)
417  msh => data%msh(1)
418 
419  do i = 1, data%size()
420  fp(i)%ptr => data%items(i)%ptr
421  end do
422  else
423  call neko_error('Empty field list')
424  end if
425 
426  nullify(dtlag)
427  nullify(tlag)
428 
429  type is(chkp_t)
430 
431  if ( .not. associated(data%u) .or. &
432  .not. associated(data%v) .or. &
433  .not. associated(data%w) .or. &
434  .not. associated(data%p) ) then
435  call neko_error('Checkpoint not initialized')
436  end if
437 
438  fp_size = 4
439 
440  if (associated(data%s)) then
441  fp_size = fp_size + 1
442  end if
443 
444  if (associated(data%abx1)) then
445  fp_size = fp_size + 6
446  end if
447 
448  if (associated(data%abs1)) then
449  fp_size = fp_size + 2
450  end if
451 
452  allocate(fp(fp_size))
453 
454  fsp_size = 0
455  if (associated(data%ulag)) then
456  fsp_size = fsp_size + 3
457  end if
458 
459  if (associated(data%slag)) then
460  fsp_size = fsp_size + 1
461  end if
462 
463  if (fsp_size .gt. 0) then
464  allocate(fsp(fsp_size))
465  fsp_cur = 1
466  end if
467 
468  dof => data%u%dof
469  msh => data%u%msh
470 
471  fp(1)%ptr => data%u
472  fp(2)%ptr => data%v
473  fp(3)%ptr => data%w
474  fp(4)%ptr => data%p
475 
476  fp_cur = 5
477  if (associated(data%s)) then
478  fp(fp_cur)%ptr => data%s
479  fp_cur = fp_cur + 1
480  end if
481 
482  if (associated(data%abx1)) then
483  fp(fp_cur)%ptr => data%abx1
484  fp(fp_cur+1)%ptr => data%abx2
485  fp(fp_cur+2)%ptr => data%aby1
486  fp(fp_cur+3)%ptr => data%aby2
487  fp(fp_cur+4)%ptr => data%abz1
488  fp(fp_cur+5)%ptr => data%abz2
489  fp_cur = fp_cur + 6
490  end if
491 
492  if (associated(data%abs1)) then
493  fp(fp_cur)%ptr => data%abs1
494  fp(fp_cur+1)%ptr => data%abs2
495  fp_cur = fp_cur + 2
496  end if
497 
498  if (associated(data%ulag)) then
499  fsp(fsp_cur)%ptr => data%ulag
500  fsp(fsp_cur+1)%ptr => data%vlag
501  fsp(fsp_cur+2)%ptr => data%wlag
502  fsp_cur = fsp_cur + 3
503  end if
504 
505  if (associated(data%slag)) then
506  fsp(fsp_cur)%ptr => data%slag
507  fsp_cur = fsp_cur + 1
508  end if
509 
510  if (associated(data%tlag)) then
511  tlag => data%tlag
512  dtlag => data%dtlag
513  end if
514 
515  class default
516  call neko_log%error('Invalid data')
517  end select
518 
519  end subroutine hdf5_file_determine_data
520 
521 #else
522 
524  subroutine hdf5_file_write(this, data, t)
525  class(hdf5_file_t), intent(inout) :: this
526  class(*), target, intent(in) :: data
527  real(kind=rp), intent(in), optional :: t
528  call neko_error('Neko needs to be built with HDF5 support')
529  end subroutine hdf5_file_write
530 
532  subroutine hdf5_file_read(this, data)
533  class(hdf5_file_t) :: this
534  class(*), target, intent(inout) :: data
535  call neko_error('Neko needs to be built with HDF5 support')
536  end subroutine hdf5_file_read
537 
538 
539 #endif
540 
541 end module hdf5_file
Defines a checkpoint.
Definition: checkpoint.f90:34
Definition: comm.F90:1
integer pe_rank
MPI rank.
Definition: comm.F90:26
type(mpi_comm) neko_comm
MPI communicator.
Definition: comm.F90:16
Defines a mapping of the degrees of freedom.
Definition: dofmap.f90:35
Stores a series fields.
Defines a field.
Definition: field.f90:34
HDF5 file format.
Definition: hdf5_file.F90:34
subroutine hdf5_file_read(this, data)
Read data in HDF5 format.
Definition: hdf5_file.F90:533
subroutine hdf5_file_write(this, data, t)
Write data in HDF5 format.
Definition: hdf5_file.F90:525
Logging routines.
Definition: log.f90:34
type(log_t), public neko_log
Global log stream.
Definition: log.f90:61
Defines a mesh.
Definition: mesh.f90:34
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Utilities.
Definition: utils.f90:35
pure integer function, public filename_suffix_pos(fname)
Find position (in the string) of a filename's suffix.
Definition: utils.f90:55
field_ptr_t, To easily obtain a pointer to a field
Definition: field.f90:80
field_list_t, To be able to group fields together
Definition: field_list.f90:13
field_series_ptr_t, To easily obtain a pointer to a field series
A generic file handler.
Interface for HDF5 files.
Definition: hdf5_file.F90:54