46 use mpi_f08,
only : mpi_info_null, mpi_allreduce, mpi_in_place, &
69 class(*),
target,
intent(in) :: data
70 real(kind=
rp),
intent(in),
optional :: t
71 type(
mesh_t),
pointer :: msh
75 real(kind=
rp),
pointer :: dtlag(:)
76 real(kind=
rp),
pointer :: tlag(:)
77 integer :: ierr, info, drank, i, j
78 integer(hid_t) :: plist_id, file_id, dset_id, grp_id, attr_id
79 integer(hid_t) :: filespace, memspace
80 integer(hid_t) :: h5t_neko_real
81 integer(hsize_t),
dimension(1) :: ddim, dcount, doffset
83 character(len=5) :: id_str
84 character(len=1024) :: fname
86 call hdf5_file_determine_data(
data, msh, dof, fp, fsp, dtlag, tlag)
88 if (.not. this%overwrite)
call this%increment_counter()
89 fname = trim(this%get_fname())
93 call hdf5_file_determine_real(h5t_neko_real)
95 call h5pcreate_f(h5p_file_access_f, plist_id, ierr)
96 info = mpi_info_null%mpi_val
97 call h5pset_fapl_mpio_f(plist_id,
neko_comm%mpi_val, info, ierr)
99 call h5fcreate_f(fname, h5f_acc_trunc_f, &
100 file_id, ierr, access_prp = plist_id)
102 call h5pcreate_f(h5p_dataset_xfer_f, plist_id, ierr)
103 call h5pset_dxpl_mpio_f(plist_id, h5fd_mpio_collective_f, ierr)
105 call h5screate_f(h5s_scalar_f, filespace, ierr)
109 call h5acreate_f(file_id,
"Time", h5t_neko_real, filespace, attr_id, &
110 ierr, h5p_default_f, h5p_default_f)
111 call h5awrite_f(attr_id, h5t_neko_real, t, ddim, ierr)
112 call h5aclose_f(attr_id, ierr)
115 if (
associated(dof))
then
116 call h5acreate_f(file_id,
"Lx", h5t_native_integer, filespace, attr_id, &
117 ierr, h5p_default_f, h5p_default_f)
118 call h5awrite_f(attr_id, h5t_native_integer, dof%Xh%lx, ddim, ierr)
119 call h5aclose_f(attr_id, ierr)
122 if (
associated(msh))
then
123 call h5gcreate_f(file_id,
"Mesh", grp_id, ierr, lcpl_id=h5p_default_f, &
124 gcpl_id=h5p_default_f, gapl_id=h5p_default_f)
126 call h5acreate_f(grp_id,
"Elements", h5t_native_integer, filespace, attr_id, &
127 ierr, h5p_default_f, h5p_default_f)
128 call h5awrite_f(attr_id, h5t_native_integer, msh%glb_nelv, ddim, ierr)
129 call h5aclose_f(attr_id, ierr)
131 call h5acreate_f(grp_id,
"Dimension", h5t_native_integer, filespace, attr_id, &
132 ierr, h5p_default_f, h5p_default_f)
133 call h5awrite_f(attr_id, h5t_native_integer, msh%gdim, ddim, ierr)
134 call h5aclose_f(attr_id, ierr)
136 call h5gclose_f(grp_id, ierr)
140 call h5sclose_f(filespace, ierr)
145 if (
associated(tlag) .and.
associated(dtlag))
then
146 call h5gcreate_f(file_id,
"Restart", grp_id, ierr, lcpl_id=h5p_default_f, &
147 gcpl_id=h5p_default_f, gapl_id=h5p_default_f)
158 call h5screate_simple_f(drank, ddim, filespace, ierr)
160 call h5dcreate_f(grp_id,
'tlag', h5t_neko_real, &
161 filespace, dset_id, ierr)
162 call h5dget_space_f(dset_id, filespace, ierr)
163 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
164 doffset, dcount, ierr)
165 call h5dwrite_f(dset_id, h5t_neko_real, tlag, &
166 ddim, ierr, xfer_prp = plist_id)
167 call h5dclose_f(dset_id, ierr)
169 call h5dcreate_f(grp_id,
'dtlag', h5t_neko_real, &
170 filespace, dset_id, ierr)
171 call h5dget_space_f(dset_id, filespace, ierr)
172 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
173 doffset, dcount, ierr)
174 call h5dwrite_f(dset_id, h5t_neko_real, dtlag, &
175 ddim, ierr, xfer_prp = plist_id)
176 call h5dclose_f(dset_id, ierr)
178 call h5sclose_f(filespace, ierr)
179 call h5gclose_f(grp_id, ierr)
187 if (
allocated(fp) .or.
allocated(fsp))
then
188 call h5gcreate_f(file_id,
"Fields", grp_id, ierr, lcpl_id=h5p_default_f, &
189 gcpl_id=h5p_default_f, gapl_id=h5p_default_f)
191 dcount(1) = int(dof%size(), 8)
192 doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
193 ddim = int(dof%size(), 8)
195 call mpi_allreduce(mpi_in_place, ddim(1), 1, &
198 call h5screate_simple_f(drank, ddim, filespace, ierr)
199 call h5screate_simple_f(drank, dcount, memspace, ierr)
202 if (
allocated(fp))
then
204 call h5dcreate_f(grp_id, fp(i)%ptr%name, h5t_neko_real, &
205 filespace, dset_id, ierr)
206 call h5dget_space_f(dset_id, filespace, ierr)
207 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
208 doffset, dcount, ierr)
209 call h5dwrite_f(dset_id, h5t_neko_real, &
210 fp(i)%ptr%x(1,1,1,1), &
211 ddim, ierr, file_space_id = filespace, &
212 mem_space_id = memspace, xfer_prp = plist_id)
213 call h5dclose_f(dset_id, ierr)
218 if (
allocated(fsp))
then
220 do j = 1, fsp(i)%ptr%size()
221 call h5dcreate_f(grp_id, fsp(i)%ptr%lf(j)%name, &
222 h5t_neko_real, filespace, dset_id, ierr)
223 call h5dget_space_f(dset_id, filespace, ierr)
224 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
225 doffset, dcount, ierr)
226 call h5dwrite_f(dset_id, h5t_neko_real, &
227 fsp(i)%ptr%lf(j)%x(1,1,1,1), &
228 ddim, ierr, file_space_id = filespace, &
229 mem_space_id = memspace, xfer_prp = plist_id)
230 call h5dclose_f(dset_id, ierr)
236 call h5gclose_f(grp_id, ierr)
237 call h5sclose_f(filespace, ierr)
238 call h5sclose_f(memspace, ierr)
241 call h5pclose_f(plist_id, ierr)
242 call h5fclose_f(file_id, ierr)
251 class(*),
target,
intent(inout) :: data
252 integer(hid_t) :: plist_id, file_id, dset_id, grp_id, attr_id
253 integer(hid_t) :: filespace, memspace
254 integer(hid_t) :: h5t_neko_real
255 integer(hsize_t),
dimension(1) :: ddim, dcount, doffset
256 integer :: i,j, ierr, info, glb_nelv, gdim, lx, drank
257 type(
mesh_t),
pointer :: msh
261 real(kind=
rp),
pointer :: dtlag(:)
262 real(kind=
rp),
pointer :: tlag(:)
264 character(len=1024) :: fname
266 fname = trim(this%get_fname())
268 call hdf5_file_determine_data(
data, msh, dof, fp, fsp, dtlag, tlag)
272 call hdf5_file_determine_real(h5t_neko_real)
274 call h5pcreate_f(h5p_file_access_f, plist_id, ierr)
275 info = mpi_info_null%mpi_val
276 call h5pset_fapl_mpio_f(plist_id,
neko_comm%mpi_val, info, ierr)
278 call h5fopen_f(fname, h5f_acc_rdonly_f, &
279 file_id, ierr, access_prp = plist_id)
281 call h5pcreate_f(h5p_dataset_xfer_f, plist_id, ierr)
282 call h5pset_dxpl_mpio_f(plist_id, h5fd_mpio_collective_f, ierr)
285 call h5aopen_name_f(file_id,
'Time', attr_id, ierr)
286 call h5aread_f(attr_id, h5t_neko_real, t, ddim, ierr)
287 call h5aclose_f(attr_id, ierr)
294 call h5aopen_name_f(file_id,
'Lx', attr_id, ierr)
295 call h5aread_f(attr_id, h5t_native_integer, lx, ddim, ierr)
296 call h5aclose_f(attr_id, ierr)
298 call h5gopen_f(file_id,
'Mesh', grp_id, ierr, gapl_id=h5p_default_f)
300 call h5aopen_name_f(grp_id,
'Elements', attr_id, ierr)
301 call h5aread_f(attr_id, h5t_native_integer, glb_nelv, ddim, ierr)
302 call h5aclose_f(attr_id, ierr)
304 call h5aopen_name_f(grp_id,
'Dimension', attr_id, ierr)
305 call h5aread_f(attr_id, h5t_native_integer, gdim, ddim, ierr)
306 call h5aclose_f(attr_id, ierr)
307 call h5gclose_f(grp_id, ierr)
310 if (
associated(tlag) .and.
associated(dtlag))
then
320 call h5gopen_f(file_id,
'Restart', grp_id, ierr, gapl_id=h5p_default_f)
321 call h5dopen_f(grp_id,
'tlag', dset_id, ierr)
322 call h5dget_space_f(dset_id, filespace, ierr)
323 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
324 doffset, dcount, ierr)
325 call h5dread_f(dset_id, h5t_neko_real, tlag, ddim, ierr, xfer_prp=plist_id)
326 call h5dclose_f(dset_id, ierr)
327 call h5sclose_f(filespace, ierr)
329 call h5dopen_f(grp_id,
'dtlag', dset_id, ierr)
330 call h5dget_space_f(dset_id, filespace, ierr)
331 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
332 doffset, dcount, ierr)
333 call h5dread_f(dset_id, h5t_neko_real, dtlag, ddim, ierr, xfer_prp=plist_id)
334 call h5dclose_f(dset_id, ierr)
335 call h5sclose_f(filespace, ierr)
337 call h5gclose_f(grp_id, ierr)
340 if (
allocated(fp) .or.
allocated(fsp))
then
341 call h5gopen_f(file_id,
'Fields', grp_id, ierr, gapl_id=h5p_default_f)
343 dcount(1) = int(dof%size(), 8)
344 doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
345 ddim = int(dof%size(), 8)
348 dcount(1) = int(dof%size(), 8)
349 doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
350 ddim = int(dof%size(), 8)
352 call mpi_allreduce(mpi_in_place, ddim(1), 1, &
355 call h5screate_simple_f(drank, dcount, memspace, ierr)
357 if (
allocated(fp))
then
359 call h5dopen_f(grp_id, fp(i)%ptr%name, dset_id, ierr)
360 call h5dget_space_f(dset_id, filespace, ierr)
361 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
362 doffset, dcount, ierr)
363 call h5dread_f(dset_id, h5t_neko_real, &
364 fp(i)%ptr%x(1,1,1,1), &
365 ddim, ierr, file_space_id = filespace, &
366 mem_space_id = memspace, xfer_prp=plist_id)
367 call h5dclose_f(dset_id, ierr)
368 call h5sclose_f(filespace, ierr)
372 if (
allocated(fsp))
then
374 do j = 1, fsp(i)%ptr%size()
375 call h5dopen_f(grp_id, fsp(i)%ptr%lf(j)%name, dset_id, ierr)
376 call h5dget_space_f(dset_id, filespace, ierr)
377 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
378 doffset, dcount, ierr)
379 call h5dread_f(dset_id, h5t_neko_real, &
380 fsp(i)%ptr%lf(j)%x(1,1,1,1), &
381 ddim, ierr, file_space_id = filespace, &
382 mem_space_id = memspace, xfer_prp=plist_id)
383 call h5dclose_f(dset_id, ierr)
384 call h5sclose_f(filespace, ierr)
388 call h5sclose_f(memspace, ierr)
389 call h5gclose_f(grp_id, ierr)
392 call h5pclose_f(plist_id, ierr)
393 call h5fclose_f(file_id, ierr)
400 subroutine hdf5_file_determine_data(data, msh, dof, fp, fsp, dtlag, tlag)
401 class(*),
target,
intent(in) :: data
402 type(
mesh_t),
pointer,
intent(inout) :: msh
403 type(
dofmap_t),
pointer,
intent(inout) :: dof
404 type(
field_ptr_t),
allocatable,
intent(inout) :: fp(:)
406 real(kind=
rp),
pointer,
intent(inout) :: dtlag(:)
407 real(kind=
rp),
pointer,
intent(inout) :: tlag(:)
408 integer :: i, j, fp_size, fp_cur, fsp_size, fsp_cur, scalar_count, ab_count
409 character(len=32) :: scalar_name
416 allocate(fp(fp_size))
424 if (data%size() .gt. 0)
then
425 allocate(fp(data%size()))
430 do i = 1, data%size()
431 fp(i)%ptr => data%items(i)%ptr
442 if ( .not.
associated(data%u) .or. &
443 .not.
associated(data%v) .or. &
444 .not.
associated(data%w) .or. &
445 .not.
associated(data%p) )
then
451 if (
allocated(data%scalar_lags%items) .and. data%scalar_lags%size() > 0)
then
452 scalar_count = data%scalar_lags%size()
453 else if (
associated(data%s))
then
459 if (scalar_count .gt. 1)
then
460 fp_size = fp_size + scalar_count
463 fp_size = fp_size + (scalar_count * 2)
464 else if (
associated(data%s))
then
466 fp_size = fp_size + 1
467 if (
associated(data%abs1))
then
468 fp_size = fp_size + 2
472 if (
associated(data%abx1))
then
473 fp_size = fp_size + 6
476 allocate(fp(fp_size))
479 if (
associated(data%ulag))
then
480 fsp_size = fsp_size + 3
483 if (scalar_count .gt. 1)
then
484 if (
allocated(data%scalar_lags%items))
then
485 fsp_size = fsp_size + data%scalar_lags%size()
487 else if (
associated(data%slag))
then
488 fsp_size = fsp_size + 1
491 if (fsp_size .gt. 0)
then
492 allocate(fsp(fsp_size))
506 if (scalar_count .gt. 1)
then
509 do i = 1, scalar_count
510 slag => data%scalar_lags%get(i)
511 fp(fp_cur)%ptr => slag%f
516 do i = 1, scalar_count
517 fp(fp_cur)%ptr => data%scalar_abx1(i)%ptr
519 fp(fp_cur)%ptr => data%scalar_abx2(i)%ptr
522 else if (
associated(data%s))
then
524 fp(fp_cur)%ptr => data%s
527 if (
associated(data%abs1))
then
528 fp(fp_cur)%ptr => data%abs1
529 fp(fp_cur+1)%ptr => data%abs2
534 if (
associated(data%abx1))
then
535 fp(fp_cur)%ptr => data%abx1
536 fp(fp_cur+1)%ptr => data%abx2
537 fp(fp_cur+2)%ptr => data%aby1
538 fp(fp_cur+3)%ptr => data%aby2
539 fp(fp_cur+4)%ptr => data%abz1
540 fp(fp_cur+5)%ptr => data%abz2
544 if (
associated(data%ulag))
then
545 fsp(fsp_cur)%ptr => data%ulag
546 fsp(fsp_cur+1)%ptr => data%vlag
547 fsp(fsp_cur+2)%ptr => data%wlag
548 fsp_cur = fsp_cur + 3
552 if (scalar_count .gt. 1)
then
553 if (
allocated(data%scalar_lags%items))
then
554 do j = 1, data%scalar_lags%size()
555 fsp(fsp_cur)%ptr => data%scalar_lags%get(j)
556 fsp_cur = fsp_cur + 1
559 else if (
associated(data%slag))
then
560 fsp(fsp_cur)%ptr => data%slag
561 fsp_cur = fsp_cur + 1
564 if (
associated(data%tlag))
then
573 end subroutine hdf5_file_determine_data
578 subroutine hdf5_file_determine_real(H5T_NEKO_REAL)
579 integer(hid_t),
intent(inout) :: h5t_neko_real
582 h5t_neko_real = h5t_native_double
584 h5t_neko_real = h5t_native_real
586 call neko_error(
"Unsupported real type")
588 end subroutine hdf5_file_determine_real
593 logical,
intent(in) :: overwrite
594 this%overwrite = overwrite
602 class(*),
target,
intent(in) :: data
603 real(kind=rp),
intent(in),
optional :: t
604 call neko_error(
'Neko needs to be built with HDF5 support')
610 class(*),
target,
intent(inout) :: data
611 call neko_error(
'Neko needs to be built with HDF5 support')
617 logical,
intent(in) :: overwrite
618 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.
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 dp
integer, parameter, public sp
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.