47 use mpi_f08,
only : mpi_info_null, mpi_allreduce, mpi_in_place, &
57 logical :: overwrite = .false.
71 class(*),
target,
intent(in) :: data
72 real(kind=
rp),
intent(in),
optional :: t
73 type(
mesh_t),
pointer :: msh
77 real(kind=
rp),
pointer :: dtlag(:)
78 real(kind=
rp),
pointer :: tlag(:)
79 integer :: ierr, info, drank, i, j
80 integer(hid_t) :: plist_id, file_id, dset_id, grp_id, attr_id
81 integer(hid_t) :: filespace, memspace
82 integer(hsize_t),
dimension(1) :: ddim, dcount, doffset
84 character(len=5) :: id_str
85 character(len=1024) :: fname
87 call hdf5_file_determine_data(
data, msh, dof, fp, fsp, dtlag, tlag)
89 if (this%overwrite)
then
90 fname = trim(this%fname)
93 write(id_str,
'(i5.5)') this%counter
94 fname = trim(this%fname(1:suffix_pos-1))//id_str//
'.h5'
98 call h5pcreate_f(h5p_file_access_f, plist_id, ierr)
99 info = mpi_info_null%mpi_val
100 call h5pset_fapl_mpio_f(plist_id,
neko_comm%mpi_val, info, ierr)
102 call h5fcreate_f(fname, h5f_acc_trunc_f, &
103 file_id, ierr, access_prp = plist_id)
105 call h5pcreate_f(h5p_dataset_xfer_f, plist_id, ierr)
106 call h5pset_dxpl_mpio_f(plist_id, h5fd_mpio_collective_f, ierr)
108 call h5screate_f(h5s_scalar_f, filespace, ierr)
112 call h5acreate_f(file_id,
"Time", h5t_native_double, filespace, attr_id, &
113 ierr, h5p_default_f, h5p_default_f)
114 call h5awrite_f(attr_id, h5t_native_double, t, ddim, ierr)
115 call h5aclose_f(attr_id, ierr)
118 if (
associated(dof))
then
119 call h5acreate_f(file_id,
"Lx", h5t_native_integer, filespace, attr_id, &
120 ierr, h5p_default_f, h5p_default_f)
121 call h5awrite_f(attr_id, h5t_native_integer, dof%Xh%lx, ddim, ierr)
122 call h5aclose_f(attr_id, ierr)
125 if (
associated(msh))
then
126 call h5gcreate_f(file_id,
"Mesh", grp_id, ierr, lcpl_id=h5p_default_f, &
127 gcpl_id=h5p_default_f, gapl_id=h5p_default_f)
129 call h5acreate_f(grp_id,
"Elements", h5t_native_integer, filespace, attr_id, &
130 ierr, h5p_default_f, h5p_default_f)
131 call h5awrite_f(attr_id, h5t_native_integer, msh%glb_nelv, ddim, ierr)
132 call h5aclose_f(attr_id, ierr)
134 call h5acreate_f(grp_id,
"Dimension", h5t_native_integer, filespace, attr_id, &
135 ierr, h5p_default_f, h5p_default_f)
136 call h5awrite_f(attr_id, h5t_native_integer, msh%gdim, ddim, ierr)
137 call h5aclose_f(attr_id, ierr)
139 call h5gclose_f(grp_id, ierr)
143 call h5sclose_f(filespace, ierr)
148 if (
associated(tlag) .and.
associated(dtlag))
then
149 call h5gcreate_f(file_id,
"Restart", grp_id, ierr, lcpl_id=h5p_default_f, &
150 gcpl_id=h5p_default_f, gapl_id=h5p_default_f)
161 call h5screate_simple_f(drank, ddim, filespace, ierr)
163 call h5dcreate_f(grp_id,
'tlag', h5t_native_double, &
164 filespace, dset_id, ierr)
165 call h5dget_space_f(dset_id, filespace, ierr)
166 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
167 doffset, dcount, ierr)
168 call h5dwrite_f(dset_id, h5t_native_double, tlag, &
169 ddim, ierr, xfer_prp = plist_id)
170 call h5dclose_f(dset_id, ierr)
172 call h5dcreate_f(grp_id,
'dtlag', h5t_native_double, &
173 filespace, dset_id, ierr)
174 call h5dget_space_f(dset_id, filespace, ierr)
175 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
176 doffset, dcount, ierr)
177 call h5dwrite_f(dset_id, h5t_native_double, dtlag, &
178 ddim, ierr, xfer_prp = plist_id)
179 call h5dclose_f(dset_id, ierr)
181 call h5sclose_f(filespace, ierr)
182 call h5gclose_f(grp_id, ierr)
190 if (
allocated(fp) .or.
allocated(fsp))
then
191 call h5gcreate_f(file_id,
"Fields", grp_id, ierr, lcpl_id=h5p_default_f, &
192 gcpl_id=h5p_default_f, gapl_id=h5p_default_f)
194 dcount(1) = int(dof%size(), 8)
195 doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
196 ddim = int(dof%size(), 8)
198 call mpi_allreduce(mpi_in_place, ddim(1), 1, &
201 call h5screate_simple_f(drank, ddim, filespace, ierr)
202 call h5screate_simple_f(drank, dcount, memspace, ierr)
205 if (
allocated(fp))
then
207 call h5dcreate_f(grp_id, fp(i)%ptr%name, h5t_native_double, &
208 filespace, dset_id, ierr)
209 call h5dget_space_f(dset_id, filespace, ierr)
210 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
211 doffset, dcount, ierr)
212 call h5dwrite_f(dset_id, h5t_native_double, &
213 fp(i)%ptr%x(1,1,1,1), &
214 ddim, ierr, file_space_id = filespace, &
215 mem_space_id = memspace, xfer_prp = plist_id)
216 call h5dclose_f(dset_id, ierr)
221 if (
allocated(fsp))
then
223 do j = 1, fsp(i)%ptr%size()
224 call h5dcreate_f(grp_id, fsp(i)%ptr%lf(j)%name, &
225 h5t_native_double, filespace, dset_id, ierr)
226 call h5dget_space_f(dset_id, filespace, ierr)
227 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
228 doffset, dcount, ierr)
229 call h5dwrite_f(dset_id, h5t_native_double, &
230 fsp(i)%ptr%lf(j)%x(1,1,1,1), &
231 ddim, ierr, file_space_id = filespace, &
232 mem_space_id = memspace, xfer_prp = plist_id)
233 call h5dclose_f(dset_id, ierr)
239 call h5gclose_f(grp_id, ierr)
240 call h5sclose_f(filespace, ierr)
241 call h5sclose_f(memspace, ierr)
244 call h5pclose_f(plist_id, ierr)
245 call h5fclose_f(file_id, ierr)
249 this%counter = this%counter + 1
256 class(*),
target,
intent(inout) :: data
257 integer(hid_t) :: plist_id, file_id, dset_id, grp_id, attr_id
258 integer(hid_t) :: filespace, memspace
259 integer(hsize_t),
dimension(1) :: ddim, dcount, doffset
260 integer :: i,j, ierr, info, glb_nelv, gdim, lx, drank
261 type(
mesh_t),
pointer :: msh
265 real(kind=
rp),
pointer :: dtlag(:)
266 real(kind=
rp),
pointer :: tlag(:)
269 call hdf5_file_determine_data(
data, msh, dof, fp, fsp, dtlag, tlag)
272 call h5pcreate_f(h5p_file_access_f, plist_id, ierr)
273 info = mpi_info_null%mpi_val
274 call h5pset_fapl_mpio_f(plist_id,
neko_comm%mpi_val, info, ierr)
276 call h5fopen_f(trim(this%fname), h5f_acc_rdonly_f, &
277 file_id, ierr, access_prp = plist_id)
279 call h5pcreate_f(h5p_dataset_xfer_f, plist_id, ierr)
280 call h5pset_dxpl_mpio_f(plist_id, h5fd_mpio_collective_f, ierr)
283 call h5aopen_name_f(file_id,
'Time', attr_id, ierr)
284 call h5aread_f(attr_id, h5t_native_double, t, ddim, ierr)
285 call h5aclose_f(attr_id, ierr)
292 call h5aopen_name_f(file_id,
'Lx', attr_id, ierr)
293 call h5aread_f(attr_id, h5t_native_integer, lx, ddim, ierr)
294 call h5aclose_f(attr_id, ierr)
296 call h5gopen_f(file_id,
'Mesh', grp_id, ierr, gapl_id=h5p_default_f)
298 call h5aopen_name_f(grp_id,
'Elements', attr_id, ierr)
299 call h5aread_f(attr_id, h5t_native_integer, glb_nelv, ddim, ierr)
300 call h5aclose_f(attr_id, ierr)
302 call h5aopen_name_f(grp_id,
'Dimension', attr_id, ierr)
303 call h5aread_f(attr_id, h5t_native_integer, gdim, ddim, ierr)
304 call h5aclose_f(attr_id, ierr)
305 call h5gclose_f(grp_id, ierr)
308 if (
associated(tlag) .and.
associated(dtlag))
then
318 call h5gopen_f(file_id,
'Restart', grp_id, ierr, gapl_id=h5p_default_f)
319 call h5dopen_f(grp_id,
'tlag', 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, tlag, ddim, ierr, xfer_prp=plist_id)
324 call h5dclose_f(dset_id, ierr)
325 call h5sclose_f(filespace, ierr)
327 call h5dopen_f(grp_id,
'dtlag', dset_id, ierr)
328 call h5dget_space_f(dset_id, filespace, ierr)
329 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
330 doffset, dcount, ierr)
331 call h5dread_f(dset_id, h5t_native_double, dtlag, ddim, ierr, xfer_prp=plist_id)
332 call h5dclose_f(dset_id, ierr)
333 call h5sclose_f(filespace, ierr)
335 call h5gclose_f(grp_id, ierr)
338 if (
allocated(fp) .or.
allocated(fsp))
then
339 call h5gopen_f(file_id,
'Fields', grp_id, ierr, gapl_id=h5p_default_f)
341 dcount(1) = int(dof%size(), 8)
342 doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
343 ddim = int(dof%size(), 8)
346 dcount(1) = int(dof%size(), 8)
347 doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
348 ddim = int(dof%size(), 8)
350 call mpi_allreduce(mpi_in_place, ddim(1), 1, &
353 call h5screate_simple_f(drank, dcount, memspace, ierr)
355 if (
allocated(fp))
then
357 call h5dopen_f(grp_id, fp(i)%ptr%name, dset_id, ierr)
358 call h5dget_space_f(dset_id, filespace, ierr)
359 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
360 doffset, dcount, ierr)
361 call h5dread_f(dset_id, h5t_native_double, &
362 fp(i)%ptr%x(1,1,1,1), &
363 ddim, ierr, file_space_id = filespace, &
364 mem_space_id = memspace, xfer_prp=plist_id)
365 call h5dclose_f(dset_id, ierr)
366 call h5sclose_f(filespace, ierr)
370 if (
allocated(fsp))
then
372 do j = 1, fsp(i)%ptr%size()
373 call h5dopen_f(grp_id, fsp(i)%ptr%lf(j)%name, dset_id, ierr)
374 call h5dget_space_f(dset_id, filespace, ierr)
375 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
376 doffset, dcount, ierr)
377 call h5dread_f(dset_id, h5t_native_double, &
378 fsp(i)%ptr%lf(j)%x(1,1,1,1), &
379 ddim, ierr, file_space_id = filespace, &
380 mem_space_id = memspace, xfer_prp=plist_id)
381 call h5dclose_f(dset_id, ierr)
382 call h5sclose_f(filespace, ierr)
386 call h5sclose_f(memspace, ierr)
387 call h5gclose_f(grp_id, ierr)
390 call h5pclose_f(plist_id, ierr)
391 call h5fclose_f(file_id, ierr)
398 subroutine hdf5_file_determine_data(data, msh, dof, fp, fsp, dtlag, tlag)
399 class(*),
target,
intent(in) :: data
400 type(
mesh_t),
pointer,
intent(inout) :: msh
401 type(
dofmap_t),
pointer,
intent(inout) :: dof
402 type(
field_ptr_t),
allocatable,
intent(inout) :: fp(:)
404 real(kind=
rp),
pointer,
intent(inout) :: dtlag(:)
405 real(kind=
rp),
pointer,
intent(inout) :: tlag(:)
406 integer :: i, j, fp_size, fp_cur, fsp_size, fsp_cur, scalar_count, ab_count
407 character(len=32) :: scalar_name
414 allocate(fp(fp_size))
422 if (data%size() .gt. 0)
then
423 allocate(fp(data%size()))
428 do i = 1, data%size()
429 fp(i)%ptr => data%items(i)%ptr
440 if ( .not.
associated(data%u) .or. &
441 .not.
associated(data%v) .or. &
442 .not.
associated(data%w) .or. &
443 .not.
associated(data%p) )
then
449 if (
allocated(data%scalar_lags%items) .and. data%scalar_lags%size() > 0)
then
450 scalar_count = data%scalar_lags%size()
451 else if (
associated(data%s))
then
457 if (scalar_count .gt. 1)
then
458 fp_size = fp_size + scalar_count
461 fp_size = fp_size + (scalar_count * 2)
462 else if (
associated(data%s))
then
464 fp_size = fp_size + 1
465 if (
associated(data%abs1))
then
466 fp_size = fp_size + 2
470 if (
associated(data%abx1))
then
471 fp_size = fp_size + 6
474 allocate(fp(fp_size))
477 if (
associated(data%ulag))
then
478 fsp_size = fsp_size + 3
481 if (scalar_count .gt. 1)
then
482 if (
allocated(data%scalar_lags%items))
then
483 fsp_size = fsp_size + data%scalar_lags%size()
485 else if (
associated(data%slag))
then
486 fsp_size = fsp_size + 1
489 if (fsp_size .gt. 0)
then
490 allocate(fsp(fsp_size))
504 if (scalar_count .gt. 1)
then
505 do i = 1, scalar_count
506 associate(slag => data%scalar_lags%get(i))
507 fp(fp_cur)%ptr => slag%f
512 do i = 1, scalar_count
513 fp(fp_cur)%ptr => data%scalar_abx1(i)%ptr
515 fp(fp_cur)%ptr => data%scalar_abx2(i)%ptr
518 else if (
associated(data%s))
then
520 fp(fp_cur)%ptr => data%s
523 if (
associated(data%abs1))
then
524 fp(fp_cur)%ptr => data%abs1
525 fp(fp_cur+1)%ptr => data%abs2
530 if (
associated(data%abx1))
then
531 fp(fp_cur)%ptr => data%abx1
532 fp(fp_cur+1)%ptr => data%abx2
533 fp(fp_cur+2)%ptr => data%aby1
534 fp(fp_cur+3)%ptr => data%aby2
535 fp(fp_cur+4)%ptr => data%abz1
536 fp(fp_cur+5)%ptr => data%abz2
540 if (
associated(data%ulag))
then
541 fsp(fsp_cur)%ptr => data%ulag
542 fsp(fsp_cur+1)%ptr => data%vlag
543 fsp(fsp_cur+2)%ptr => data%wlag
544 fsp_cur = fsp_cur + 3
548 if (scalar_count .gt. 1)
then
549 if (
allocated(data%scalar_lags%items))
then
550 do j = 1, data%scalar_lags%size()
551 fsp(fsp_cur)%ptr => data%scalar_lags%get(j)
552 fsp_cur = fsp_cur + 1
555 else if (
associated(data%slag))
then
556 fsp(fsp_cur)%ptr => data%slag
557 fsp_cur = fsp_cur + 1
560 if (
associated(data%tlag))
then
569 end subroutine hdf5_file_determine_data
574 logical,
intent(in) :: overwrite
575 this%overwrite = overwrite
583 class(*),
target,
intent(in) :: data
584 real(kind=rp),
intent(in),
optional :: t
585 call neko_error(
'Neko needs to be built with HDF5 support')
591 class(*),
target,
intent(inout) :: data
592 call neko_error(
'Neko needs to be built with HDF5 support')
598 logical,
intent(in) :: overwrite
599 call neko_error(
'Neko needs to be built with HDF5 support')
integer, public pe_rank
MPI rank.
type(mpi_comm), public neko_comm
MPI communicator.
Defines a mapping of the degrees of freedom.
Defines a registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
Contains the field_serties_t type.
subroutine hdf5_file_read(this, data)
Read data in HDF5 format.
subroutine hdf5_file_write(this, data, t)
Write data in HDF5 format.
subroutine hdf5_file_set_overwrite(this, overwrite)
Set the overwrite flag for HDF5 files.
type(log_t), public neko_log
Global log stream.
integer, parameter, public rp
Global precision used in computations.
pure integer function, public filename_suffix_pos(fname)
Find position (in the string) of a filename's suffix.
field_ptr_t, To easily obtain a pointer to a field
field_list_t, To be able to group fields together
A wrapper for a pointer to a field_series_t.
Stores a series (sequence) of fields, logically connected to a base field, and arranged according to ...
Interface for HDF5 files.