50 use mpi_f08,
only : mpi_info_null, mpi_allreduce, mpi_allgather, &
51 mpi_in_place, mpi_integer, mpi_sum, mpi_max, mpi_comm_size, mpi_exscan, &
52 mpi_barrier, mpi_integer8, mpi_scan
64 integer(hid_t) :: file_id = -1_hid_t
65 integer(hid_t) :: active_group_id = -1_hid_t
66 integer(hid_t) :: plist_id = -1_hid_t
68 character(len=1) :: mode
69 integer :: precision = -1
86 procedure, pass(this) :: write_int_attribute => &
104 logical,
intent(in) :: overwrite
105 this%overwrite = overwrite
111 character(len=1024) :: base_fname
112 character(len=1024) :: fname
113 character(len=1024) :: path, name, suffix
115 fname = trim(this%get_base_fname())
123 base_fname = trim(fname)
130 integer,
intent(in) :: precision
131 this%precision = precision
144 class(*),
target,
intent(in) :: data
145 real(kind=
rp),
intent(in),
optional :: t
146 type(
mesh_t),
pointer :: msh
150 real(kind=
rp),
pointer :: dtlag(:)
151 real(kind=
rp),
pointer :: tlag(:)
152 integer :: ierr, info, drank, i, j
153 integer(hid_t) :: plist_id, file_id, dset_id, grp_id, attr_id
154 integer(hid_t) :: filespace, memspace
155 integer(hid_t) :: H5T_NEKO_REAL
156 integer(hsize_t),
dimension(1) :: ddim, dcount, doffset
157 integer :: suffix_pos
158 character(len=5) :: id_str
159 character(len=1024) :: fname
161 call hdf5_file_determine_data(
data, msh, dof, fp, fsp, dtlag, tlag)
163 if (.not. this%overwrite)
call this%increment_counter()
164 fname = trim(this%get_fname())
168 call hdf5_file_determine_real(h5t_neko_real)
170 call h5pcreate_f(h5p_file_access_f, plist_id, ierr)
171 info = mpi_info_null%mpi_val
172 call h5pset_fapl_mpio_f(plist_id,
neko_comm%mpi_val, info, ierr)
174 call h5fcreate_f(fname, h5f_acc_trunc_f, &
175 file_id, ierr, access_prp = plist_id)
177 call h5pcreate_f(h5p_dataset_xfer_f, plist_id, ierr)
178 call h5pset_dxpl_mpio_f(plist_id, h5fd_mpio_collective_f, ierr)
180 call h5screate_f(h5s_scalar_f, filespace, ierr)
184 call h5acreate_f(file_id,
"Time", h5t_neko_real, filespace, attr_id, &
185 ierr, h5p_default_f, h5p_default_f)
186 call h5awrite_f(attr_id, h5t_neko_real, t, ddim, ierr)
187 call h5aclose_f(attr_id, ierr)
190 if (
associated(dof))
then
191 call h5acreate_f(file_id,
"Lx", h5t_native_integer, filespace, attr_id, &
192 ierr, h5p_default_f, h5p_default_f)
193 call h5awrite_f(attr_id, h5t_native_integer, dof%Xh%lx, ddim, ierr)
194 call h5aclose_f(attr_id, ierr)
197 if (
associated(msh))
then
198 call h5gcreate_f(file_id,
"Mesh", grp_id, ierr, &
199 lcpl_id = h5p_default_f, gcpl_id = h5p_default_f, &
200 gapl_id = h5p_default_f)
202 call h5acreate_f(grp_id,
"Elements", h5t_native_integer, filespace, &
203 attr_id, ierr, h5p_default_f, h5p_default_f)
204 call h5awrite_f(attr_id, h5t_native_integer, msh%glb_nelv, ddim, ierr)
205 call h5aclose_f(attr_id, ierr)
207 call h5acreate_f(grp_id,
"Dimension", h5t_native_integer, filespace, &
208 attr_id, ierr, h5p_default_f, h5p_default_f)
209 call h5awrite_f(attr_id, h5t_native_integer, msh%gdim, ddim, ierr)
210 call h5aclose_f(attr_id, ierr)
212 call h5gclose_f(grp_id, ierr)
216 call h5sclose_f(filespace, ierr)
221 if (
associated(tlag) .and.
associated(dtlag))
then
222 call h5gcreate_f(file_id,
"Restart", grp_id, ierr, &
223 lcpl_id = h5p_default_f, gcpl_id = h5p_default_f, &
224 gapl_id = h5p_default_f)
235 call h5screate_simple_f(drank, ddim, filespace, ierr)
237 call h5dcreate_f(grp_id,
'tlag', h5t_neko_real, &
238 filespace, dset_id, ierr)
239 call h5dget_space_f(dset_id, filespace, ierr)
240 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
241 doffset, dcount, ierr)
242 call h5dwrite_f(dset_id, h5t_neko_real, tlag, &
243 ddim, ierr, xfer_prp = plist_id)
244 call h5dclose_f(dset_id, ierr)
246 call h5dcreate_f(grp_id,
'dtlag', h5t_neko_real, &
247 filespace, dset_id, ierr)
248 call h5dget_space_f(dset_id, filespace, ierr)
249 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
250 doffset, dcount, ierr)
251 call h5dwrite_f(dset_id, h5t_neko_real, dtlag, &
252 ddim, ierr, xfer_prp = plist_id)
253 call h5dclose_f(dset_id, ierr)
255 call h5sclose_f(filespace, ierr)
256 call h5gclose_f(grp_id, ierr)
264 if (
allocated(fp) .or.
allocated(fsp))
then
265 call h5gcreate_f(file_id,
"Fields", grp_id, ierr, &
266 lcpl_id = h5p_default_f, gcpl_id = h5p_default_f, &
267 gapl_id = h5p_default_f)
269 dcount(1) = int(dof%size(), 8)
270 doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
271 ddim = int(dof%size(), 8)
273 call mpi_allreduce(mpi_in_place, ddim(1), 1, &
276 call h5screate_simple_f(drank, ddim, filespace, ierr)
277 call h5screate_simple_f(drank, dcount, memspace, ierr)
280 if (
allocated(fp))
then
282 call h5dcreate_f(grp_id, fp(i)%ptr%name, h5t_neko_real, &
283 filespace, dset_id, ierr)
284 call h5dget_space_f(dset_id, filespace, ierr)
285 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
286 doffset, dcount, ierr)
287 call h5dwrite_f(dset_id, h5t_neko_real, &
288 fp(i)%ptr%x(1,1,1,1), &
289 ddim, ierr, file_space_id = filespace, &
290 mem_space_id = memspace, xfer_prp = plist_id)
291 call h5dclose_f(dset_id, ierr)
296 if (
allocated(fsp))
then
298 do j = 1, fsp(i)%ptr%size()
299 call h5dcreate_f(grp_id, fsp(i)%ptr%lf(j)%name, &
300 h5t_neko_real, filespace, dset_id, ierr)
301 call h5dget_space_f(dset_id, filespace, ierr)
302 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
303 doffset, dcount, ierr)
304 call h5dwrite_f(dset_id, h5t_neko_real, &
305 fsp(i)%ptr%lf(j)%x(1,1,1,1), &
306 ddim, ierr, file_space_id = filespace, &
307 mem_space_id = memspace, xfer_prp = plist_id)
308 call h5dclose_f(dset_id, ierr)
314 call h5gclose_f(grp_id, ierr)
315 call h5sclose_f(filespace, ierr)
316 call h5sclose_f(memspace, ierr)
319 call h5pclose_f(plist_id, ierr)
320 call h5fclose_f(file_id, ierr)
329 class(*),
target,
intent(inout) :: data
330 integer(hid_t) :: plist_id, file_id, dset_id, grp_id, attr_id
331 integer(hid_t) :: filespace, memspace
332 integer(hid_t) :: H5T_NEKO_REAL
333 integer(hsize_t),
dimension(1) :: ddim, dcount, doffset
334 integer :: i,j, ierr, info, glb_nelv, gdim, lx, drank
335 type(
mesh_t),
pointer :: msh
339 real(kind=
rp),
pointer :: dtlag(:)
340 real(kind=
rp),
pointer :: tlag(:)
342 character(len=1024) :: fname
344 fname = trim(this%get_fname())
346 call hdf5_file_determine_data(
data, msh, dof, fp, fsp, dtlag, tlag)
350 call hdf5_file_determine_real(h5t_neko_real)
352 call h5pcreate_f(h5p_file_access_f, plist_id, ierr)
353 info = mpi_info_null%mpi_val
354 call h5pset_fapl_mpio_f(plist_id,
neko_comm%mpi_val, info, ierr)
356 call h5fopen_f(fname, h5f_acc_rdonly_f, &
357 file_id, ierr, access_prp = plist_id)
359 call h5pcreate_f(h5p_dataset_xfer_f, plist_id, ierr)
360 call h5pset_dxpl_mpio_f(plist_id, h5fd_mpio_collective_f, ierr)
363 call h5aopen_name_f(file_id,
'Time', attr_id, ierr)
364 call h5aread_f(attr_id, h5t_neko_real, t, ddim, ierr)
365 call h5aclose_f(attr_id, ierr)
372 call h5aopen_name_f(file_id,
'Lx', attr_id, ierr)
373 call h5aread_f(attr_id, h5t_native_integer, lx, ddim, ierr)
374 call h5aclose_f(attr_id, ierr)
376 call h5gopen_f(file_id,
'Mesh', grp_id, ierr, gapl_id = h5p_default_f)
378 call h5aopen_name_f(grp_id,
'Elements', attr_id, ierr)
379 call h5aread_f(attr_id, h5t_native_integer, glb_nelv, ddim, ierr)
380 call h5aclose_f(attr_id, ierr)
382 call h5aopen_name_f(grp_id,
'Dimension', attr_id, ierr)
383 call h5aread_f(attr_id, h5t_native_integer, gdim, ddim, ierr)
384 call h5aclose_f(attr_id, ierr)
385 call h5gclose_f(grp_id, ierr)
388 if (
associated(tlag) .and.
associated(dtlag))
then
398 call h5gopen_f(file_id,
'Restart', grp_id, ierr, &
399 gapl_id = h5p_default_f)
400 call h5dopen_f(grp_id,
'tlag', dset_id, ierr)
401 call h5dget_space_f(dset_id, filespace, ierr)
402 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
403 doffset, dcount, ierr)
404 call h5dread_f(dset_id, h5t_neko_real, tlag, ddim, ierr, &
406 call h5dclose_f(dset_id, ierr)
407 call h5sclose_f(filespace, ierr)
409 call h5dopen_f(grp_id,
'dtlag', dset_id, ierr)
410 call h5dget_space_f(dset_id, filespace, ierr)
411 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
412 doffset, dcount, ierr)
413 call h5dread_f(dset_id, h5t_neko_real, dtlag, ddim, ierr, &
415 call h5dclose_f(dset_id, ierr)
416 call h5sclose_f(filespace, ierr)
418 call h5gclose_f(grp_id, ierr)
421 if (
allocated(fp) .or.
allocated(fsp))
then
422 call h5gopen_f(file_id,
'Fields', grp_id, ierr, gapl_id = h5p_default_f)
424 dcount(1) = int(dof%size(), 8)
425 doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
426 ddim = int(dof%size(), 8)
429 dcount(1) = int(dof%size(), 8)
430 doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
431 ddim = int(dof%size(), 8)
433 call mpi_allreduce(mpi_in_place, ddim(1), 1, &
436 call h5screate_simple_f(drank, dcount, memspace, ierr)
438 if (
allocated(fp))
then
440 call h5dopen_f(grp_id, fp(i)%ptr%name, dset_id, ierr)
441 call h5dget_space_f(dset_id, filespace, ierr)
442 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
443 doffset, dcount, ierr)
444 call h5dread_f(dset_id, h5t_neko_real, &
445 fp(i)%ptr%x(1,1,1,1), &
446 ddim, ierr, file_space_id = filespace, &
447 mem_space_id = memspace, xfer_prp = plist_id)
448 call h5dclose_f(dset_id, ierr)
449 call h5sclose_f(filespace, ierr)
453 if (
allocated(fsp))
then
455 do j = 1, fsp(i)%ptr%size()
456 call h5dopen_f(grp_id, fsp(i)%ptr%lf(j)%name, dset_id, ierr)
457 call h5dget_space_f(dset_id, filespace, ierr)
458 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
459 doffset, dcount, ierr)
460 call h5dread_f(dset_id, h5t_neko_real, &
461 fsp(i)%ptr%lf(j)%x(1,1,1,1), &
462 ddim, ierr, file_space_id = filespace, &
463 mem_space_id = memspace, xfer_prp = plist_id)
464 call h5dclose_f(dset_id, ierr)
465 call h5sclose_f(filespace, ierr)
469 call h5sclose_f(memspace, ierr)
470 call h5gclose_f(grp_id, ierr)
473 call h5pclose_f(plist_id, ierr)
474 call h5fclose_f(file_id, ierr)
481 subroutine hdf5_file_determine_data(data, msh, dof, fp, fsp, dtlag, tlag)
482 class(*),
target,
intent(in) :: data
483 type(
mesh_t),
pointer,
intent(inout) :: msh
484 type(
dofmap_t),
pointer,
intent(inout) :: dof
485 type(
field_ptr_t),
allocatable,
intent(inout) :: fp(:)
487 real(kind=
rp),
pointer,
intent(inout) :: dtlag(:)
488 real(kind=
rp),
pointer,
intent(inout) :: tlag(:)
489 integer :: i, j, fp_size, fp_cur, fsp_size, fsp_cur, scalar_count, ab_count
490 character(len=32) :: scalar_name
497 allocate(fp(fp_size))
505 if (data%size() .gt. 0)
then
506 allocate(fp(data%size()))
511 do i = 1, data%size()
512 fp(i)%ptr => data%items(i)%ptr
523 if ( .not.
associated(data%u) .or. &
524 .not.
associated(data%v) .or. &
525 .not.
associated(data%w) .or. &
526 .not.
associated(data%p) )
then
532 if (
allocated(data%scalar_lags%items) .and. &
533 data%scalar_lags%size() > 0)
then
534 scalar_count = data%scalar_lags%size()
535 else if (
associated(data%s))
then
541 if (scalar_count .gt. 1)
then
542 fp_size = fp_size + scalar_count
545 fp_size = fp_size + (scalar_count * 2)
546 else if (
associated(data%s))
then
548 fp_size = fp_size + 1
549 if (
associated(data%abs1))
then
550 fp_size = fp_size + 2
554 if (
associated(data%abx1))
then
555 fp_size = fp_size + 6
558 allocate(fp(fp_size))
561 if (
associated(data%ulag))
then
562 fsp_size = fsp_size + 3
565 if (scalar_count .gt. 1)
then
566 if (
allocated(data%scalar_lags%items))
then
567 fsp_size = fsp_size + data%scalar_lags%size()
569 else if (
associated(data%slag))
then
570 fsp_size = fsp_size + 1
573 if (fsp_size .gt. 0)
then
574 allocate(fsp(fsp_size))
588 if (scalar_count .gt. 1)
then
591 do i = 1, scalar_count
592 slag => data%scalar_lags%get(i)
593 fp(fp_cur)%ptr => slag%f
598 do i = 1, scalar_count
599 fp(fp_cur)%ptr => data%scalar_abx1(i)%ptr
601 fp(fp_cur)%ptr => data%scalar_abx2(i)%ptr
604 else if (
associated(data%s))
then
606 fp(fp_cur)%ptr => data%s
609 if (
associated(data%abs1))
then
610 fp(fp_cur)%ptr => data%abs1
611 fp(fp_cur+1)%ptr => data%abs2
616 if (
associated(data%abx1))
then
617 fp(fp_cur)%ptr => data%abx1
618 fp(fp_cur+1)%ptr => data%abx2
619 fp(fp_cur+2)%ptr => data%aby1
620 fp(fp_cur+3)%ptr => data%aby2
621 fp(fp_cur+4)%ptr => data%abz1
622 fp(fp_cur+5)%ptr => data%abz2
626 if (
associated(data%ulag))
then
627 fsp(fsp_cur)%ptr => data%ulag
628 fsp(fsp_cur+1)%ptr => data%vlag
629 fsp(fsp_cur+2)%ptr => data%wlag
630 fsp_cur = fsp_cur + 3
634 if (scalar_count .gt. 1)
then
635 if (
allocated(data%scalar_lags%items))
then
636 do j = 1, data%scalar_lags%size()
637 fsp(fsp_cur)%ptr => data%scalar_lags%get(j)
638 fsp_cur = fsp_cur + 1
641 else if (
associated(data%slag))
then
642 fsp(fsp_cur)%ptr => data%slag
643 fsp_cur = fsp_cur + 1
646 if (
associated(data%tlag))
then
655 end subroutine hdf5_file_determine_data
660 subroutine hdf5_file_determine_real(H5T_NEKO_REAL)
661 integer(hid_t),
intent(inout) :: H5T_NEKO_REAL
664 h5t_neko_real = h5t_native_double
666 h5t_neko_real = h5t_native_real
668 call neko_error(
"Unsupported real type")
670 end subroutine hdf5_file_determine_real
679 character(len=1),
intent(in) :: mode
680 integer :: ierr, mpi_info, mpi_comm, i, n_fields, counter
681 logical :: file_exists
682 character(len=1024) :: fname
683 character(len=LOG_SIZE) :: log_buf
689 if (this%precision .gt. rp)
then
691 call neko_warning(
'Requested precision is higher than working precision')
692 else if (this%precision .eq. -1)
then
697 counter = this%get_counter() - this%get_start_counter()
700 mpi_info = mpi_info_null%mpi_val
701 mpi_comm = neko_comm%mpi_val
703 call h5pcreate_f(h5p_file_access_f, this%plist_id, ierr)
704 call h5pset_fapl_mpio_f(this%plist_id, mpi_comm, mpi_info, ierr)
707 inquire(
file = fname, exist = file_exists)
708 if (file_exists)
then
709 call h5fopen_f(fname, h5f_acc_rdwr_f, this%file_id, ierr, &
710 access_prp = this%plist_id)
712 call h5fcreate_f(fname, h5f_acc_trunc_f, &
713 this%file_id, ierr, access_prp = this%plist_id)
717 call this%set_active_group()
719 write (log_buf, *)
"Opened HDF5 file: ", trim(fname),
" with counter: ", &
721 call neko_log%message(log_buf, lvl = neko_log_debug)
730 if (this%active_group_id .ne. -1_hid_t .and. &
731 this%active_group_id .ne. this%file_id)
then
732 call h5gclose_f(this%active_group_id, ierr)
734 this%active_group_id = -1_hid_t
736 call h5pclose_f(this%plist_id, ierr)
737 this%plist_id = -1_hid_t
738 call h5fclose_f(this%file_id, ierr)
739 this%file_id = -1_hid_t
742 call neko_log%message(
"Closed HDF5 file: " // trim(this%get_fname()), &
743 lvl = neko_log_debug)
753 character(len=*),
intent(in),
optional :: group_name_path
754 character(len=1000),
allocatable :: group_name(:)
756 integer(hid_t) :: current_id, group_id
757 integer :: ierr, i, j, num_groups, name_len, group_loc
758 logical :: group_exists
762 if (this%active_group_id .ne. -1_hid_t .and. this%active_group_id .ne. &
764 call h5gclose_f(this%active_group_id, ierr)
766 this%active_group_id = -1_hid_t
769 current_id = this%file_id
771 if (.not.
present(group_name_path))
then
772 this%active_group_id = current_id
777 name_len = len(trim(group_name_path))
781 if (group_name_path .eq.
"/")
then
782 num_groups = num_groups + 1
787 allocate(group_name(num_groups))
791 if (group_name_path .eq.
"/")
then
792 group_name(group_loc) = group_name_path(j:i-1)
793 group_loc = group_loc + 1
797 if (j .ne. name_len)
then
798 group_name(group_loc) = group_name_path(j:name_len)
803 call h5lexists_f(current_id, trim(group_name(i)), group_exists, ierr)
806 if (group_exists)
then
807 call h5gopen_f(current_id, trim(group_name(i)), group_id, ierr)
809 if (this%mode ==
"r")
then
810 call neko_error(
"Group " // trim(group_name(i)) // &
813 call h5gcreate_f(current_id, trim(group_name(i)), group_id, ierr)
818 call h5gclose_f(current_id, ierr)
821 current_id = group_id
824 this%active_group_id = current_id
830 class(*),
intent(inout) :: data
832 select type (d => data)
834 call this%write_vector(d)
836 call this%write_matrix(d)
838 call this%write_field(d)
840 call neko_error(
"write_dataset not implemented for this data type")
846 character(len=*),
intent(in) :: data_name
847 class(*),
intent(inout) :: data
848 character(len=*),
intent(in),
optional :: strategy
850 select type (d => data)
852 call this%read_vector(data_name, d, strategy)
854 call this%read_matrix(data_name, d, strategy)
856 call neko_error(
"Reading a field_t is not supported yet")
858 call neko_error(
"read_dataset not implemented for this data type")
864 character(len=*),
intent(in) :: data_name
865 class(*),
intent(inout) :: data
867 select type (d => data)
869 call this%write_int_attribute(data_name, d)
870 type is (
real(kind=rp))
871 call this%write_rp_attribute(data_name, d)
873 call neko_error(
"write_attribute not implemented for this data type")
879 character(len=*),
intent(in) :: data_name
880 class(*),
intent(inout) :: data
881 logical,
intent(inout) :: exist
883 select type (d => data)
885 call this%read_int_attribute(data_name, d, exist)
886 type is (
real(kind=rp))
887 call this%read_rp_attribute(data_name, d, exist)
889 call neko_error(
"read_attribute not implemented for this data type")
896 type(vector_t),
intent(inout) :: vec
897 integer :: ierr, counts, offset, total_count, dset_rank, max_count
898 integer(hsize_t) :: append_offset
899 integer(hid_t) :: precision_hdf
900 integer(hid_t) :: xf_id, filespace, dset_id, memspace, dcpl_id
901 integer(hsize_t),
dimension(1) :: dcount, doffset
902 integer(hsize_t),
dimension(1) :: ddims, ddims_max, chunkdims
903 integer(hsize_t),
dimension(1) :: tempddims, tempmaxddims
904 logical :: dset_exists
905 real(kind=sp),
allocatable :: write_buffer_sp(:)
906 real(kind=dp),
allocatable :: write_buffer_dp(:)
912 append_offset = 0_hsize_t
916 call mpi_scan(counts, offset, 1, mpi_integer, &
917 mpi_sum, neko_comm, ierr)
918 offset = offset - counts
919 call mpi_allreduce(counts, total_count, 1, mpi_integer, &
920 mpi_sum, neko_comm, ierr)
921 call mpi_allreduce(counts, max_count, 1, mpi_integer, &
922 mpi_max, neko_comm, ierr)
927 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
928 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
929 precision_hdf = h5kind_to_type(this%precision, h5_real_kind)
935 ddims = [int(total_count, hsize_t)]
938 chunkdims = [
max(int(max_count, hsize_t), 1_hsize_t)]
939 ddims_max = [h5s_unlimited_f]
940 call h5lexists_f(this%active_group_id, trim(vec%name), dset_exists, ierr)
941 if (dset_exists)
then
942 if (this%overwrite)
then
944 call h5dopen_f(this%active_group_id, trim(vec%name), dset_id, ierr)
947 call h5dopen_f(this%active_group_id, trim(vec%name), dset_id, ierr)
949 call h5dget_space_f(dset_id, filespace, ierr)
951 call h5sget_simple_extent_dims_f(filespace, tempddims, tempmaxddims, &
954 call h5sclose_f(filespace, ierr)
956 ddims(1) = ddims(1) + tempddims(1)
957 append_offset = tempddims(1)
959 call h5dset_extent_f(dset_id, ddims, ierr)
963 call h5screate_simple_f(dset_rank, ddims, filespace, ierr, ddims_max)
965 call h5pcreate_f(h5p_dataset_create_f, dcpl_id, ierr)
966 call h5pset_chunk_f(dcpl_id, dset_rank, chunkdims, ierr)
968 call h5dcreate_f(this%active_group_id, trim(vec%name), precision_hdf, &
969 filespace, dset_id, ierr, dcpl_id = dcpl_id)
971 call h5sclose_f(filespace, ierr)
972 call h5pclose_f(dcpl_id, ierr)
978 dcount = [int(counts, hsize_t)]
981 doffset = [int(offset, hsize_t) + append_offset]
983 call h5dget_space_f(dset_id, filespace, ierr)
985 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, doffset, dcount, &
988 call h5screate_simple_f(dset_rank, dcount, memspace, ierr)
994 if (this%precision == sp)
then
995 allocate(write_buffer_sp(vec%size()))
996 if (vec%size() > 0) write_buffer_sp =
real(vec%x, kind=sp)
998 call h5dwrite_f(dset_id, precision_hdf, write_buffer_sp, dcount, ierr, &
999 file_space_id = filespace, mem_space_id = memspace, &
1001 deallocate(write_buffer_sp)
1002 else if (this%precision == dp)
then
1003 allocate(write_buffer_dp(vec%size()))
1004 if (vec%size() > 0) write_buffer_dp =
real(vec%x, kind=dp)
1006 call h5dwrite_f(dset_id, precision_hdf, write_buffer_dp, dcount, ierr, &
1007 file_space_id = filespace, mem_space_id = memspace, &
1009 deallocate(write_buffer_dp)
1011 call neko_error(
"Unsupported precision")
1017 call h5pclose_f(xf_id, ierr)
1018 call h5sclose_f(memspace, ierr)
1019 call h5sclose_f(filespace, ierr)
1020 call h5dclose_f(dset_id, ierr)
1026 type(matrix_t),
intent(inout) :: mat
1027 integer :: ierr, counts, offset, total_count, dset_rank, strides, max_count
1028 integer(hsize_t) :: append_offset
1029 integer(hid_t) :: precision_hdf
1030 integer(hid_t) :: xf_id, filespace, dset_id, memspace, dcpl_id
1031 integer(hsize_t),
dimension(2) :: dcount, doffset
1032 integer(hsize_t),
dimension(2) :: ddims, ddims_max, chunkdims
1033 integer(hsize_t),
dimension(2) :: tempddims, tempmaxddims
1034 logical :: dset_exists
1035 real(kind=sp),
allocatable :: write_buffer_sp(:,:)
1036 real(kind=dp),
allocatable :: write_buffer_dp(:,:)
1041 strides = mat%get_nrows()
1042 counts = mat%get_ncols()
1043 append_offset = 0_hsize_t
1047 call mpi_scan(counts, offset, 1, mpi_integer, &
1048 mpi_sum, neko_comm, ierr)
1049 offset = offset - counts
1050 call mpi_allreduce(counts, total_count, 1, mpi_integer, &
1051 mpi_sum, neko_comm, ierr)
1052 call mpi_allreduce(counts, max_count, 1, mpi_integer, &
1053 mpi_max, neko_comm, ierr)
1058 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
1059 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
1060 precision_hdf = h5kind_to_type(this%precision, h5_real_kind)
1067 ddims = [int(strides, hsize_t), int(total_count, hsize_t)]
1068 chunkdims = [int(strides, hsize_t),
max(int(max_count, hsize_t), 1_hsize_t)]
1069 ddims_max = [int(strides, hsize_t), h5s_unlimited_f]
1070 call h5lexists_f(this%active_group_id, trim(mat%name), dset_exists, ierr)
1071 if (dset_exists)
then
1072 if (this%overwrite)
then
1074 if (pe_rank .eq. 0)
then
1075 call neko_warning(
"Dataset " // trim(mat%name) // &
1076 " already exists and wil be overwritten")
1079 call h5dopen_f(this%active_group_id, trim(mat%name), dset_id, ierr)
1081 call h5dopen_f(this%active_group_id, trim(mat%name), dset_id, ierr)
1082 call h5dget_space_f(dset_id, filespace, ierr)
1083 call h5sget_simple_extent_dims_f(filespace, tempddims, tempmaxddims, &
1085 call h5sclose_f(filespace, ierr)
1086 ddims(2) = ddims(2) + tempddims(2)
1087 append_offset = tempddims(2)
1088 call h5dset_extent_f(dset_id, ddims, ierr)
1092 call h5screate_simple_f(dset_rank, ddims, filespace, ierr, ddims_max)
1093 call h5pcreate_f(h5p_dataset_create_f, dcpl_id, ierr)
1094 call h5pset_chunk_f(dcpl_id, dset_rank, chunkdims, ierr)
1096 call h5dcreate_f(this%active_group_id, trim(mat%name), precision_hdf, &
1097 filespace, dset_id, ierr, dcpl_id = dcpl_id)
1098 call h5sclose_f(filespace, ierr)
1099 call h5pclose_f(dcpl_id, ierr)
1106 dcount = [int(strides, hsize_t), int(counts, hsize_t)]
1108 doffset = [0_hsize_t, int(offset, hsize_t) + append_offset]
1110 call h5dget_space_f(dset_id, filespace, ierr)
1112 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, doffset, dcount, &
1115 call h5screate_simple_f(dset_rank, dcount, memspace, ierr)
1120 if (this%precision == sp)
then
1121 allocate(write_buffer_sp(mat%get_nrows(), mat%get_ncols()))
1122 if (mat%size() > 0) write_buffer_sp =
real(mat%x, kind=sp)
1124 call h5dwrite_f(dset_id, precision_hdf, write_buffer_sp, dcount, ierr, &
1125 file_space_id = filespace, mem_space_id = memspace, &
1127 deallocate(write_buffer_sp)
1128 else if (this%precision == dp)
then
1129 allocate(write_buffer_dp(mat%get_nrows(), mat%get_ncols()))
1130 if (mat%size() > 0) write_buffer_dp =
real(mat%x, kind=dp)
1132 call h5dwrite_f(dset_id, precision_hdf, write_buffer_dp, dcount, ierr, &
1133 file_space_id = filespace, mem_space_id = memspace, &
1135 deallocate(write_buffer_dp)
1137 call neko_error(
"Unsupported precision")
1143 call h5pclose_f(xf_id, ierr)
1144 call h5sclose_f(memspace, ierr)
1145 call h5sclose_f(filespace, ierr)
1146 call h5dclose_f(dset_id, ierr)
1152 type(field_t),
intent(inout) :: field
1153 integer :: ierr, counts, offset, total_count, dset_rank, max_count
1154 integer :: stride_ax_1, stride_ax_2, stride_ax_3
1155 integer(hsize_t) :: append_offset
1156 integer(hid_t) :: precision_hdf
1157 integer(hid_t) :: xf_id, filespace, dset_id, memspace, dcpl_id
1158 integer(hsize_t),
dimension(4) :: dcount, doffset
1159 integer(hsize_t),
dimension(4) :: ddims, ddims_max, chunkdims
1160 integer(hsize_t),
dimension(4) :: tempddims, tempmaxddims
1161 logical :: dset_exists
1162 real(kind=sp),
allocatable :: write_buffer_sp(:,:,:,:)
1163 real(kind=dp),
allocatable :: write_buffer_dp(:,:,:,:)
1168 stride_ax_1 =
field%Xh%lx
1169 stride_ax_2 =
field%Xh%ly
1170 stride_ax_3 =
field%Xh%lz
1171 counts =
field%msh%nelv
1172 append_offset = 0_hsize_t
1173 total_count =
field%msh%glb_nelv
1175 offset =
field%msh%offset_el
1176 call mpi_allreduce(counts, max_count, 1, mpi_integer, &
1177 mpi_max, neko_comm, ierr)
1182 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
1183 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
1184 precision_hdf = h5kind_to_type(this%precision, h5_real_kind)
1190 ddims = [int(stride_ax_1, hsize_t), &
1191 int(stride_ax_2, hsize_t), &
1192 int(stride_ax_3, hsize_t), &
1193 int(total_count, hsize_t)]
1194 chunkdims = [int(stride_ax_1, hsize_t), &
1195 int(stride_ax_2, hsize_t), &
1196 int(stride_ax_3, hsize_t), &
1197 max(int(max_count, hsize_t), 1_hsize_t)]
1198 ddims_max = [int(stride_ax_1, hsize_t), &
1199 int(stride_ax_2, hsize_t), &
1200 int(stride_ax_3, hsize_t), &
1202 call h5lexists_f(this%active_group_id, trim(
field%name), dset_exists, ierr)
1203 if (dset_exists)
then
1204 if (this%overwrite)
then
1206 if (pe_rank .eq. 0)
then
1207 call neko_warning(
"Overwriting dataset: " // trim(
field%name))
1209 call h5dopen_f(this%active_group_id, trim(
field%name), dset_id, ierr)
1211 call h5dopen_f(this%active_group_id, trim(
field%name), dset_id, ierr)
1212 call h5dget_space_f(dset_id, filespace, ierr)
1213 call h5sget_simple_extent_dims_f(filespace, tempddims, tempmaxddims, &
1215 call h5sclose_f(filespace, ierr)
1216 ddims(4) = ddims(4) + tempddims(4)
1217 append_offset = tempddims(4)
1218 call h5dset_extent_f(dset_id, ddims, ierr)
1222 call h5screate_simple_f(dset_rank, ddims, filespace, ierr, ddims_max)
1223 call h5pcreate_f(h5p_dataset_create_f, dcpl_id, ierr)
1224 call h5pset_chunk_f(dcpl_id, dset_rank, chunkdims, ierr)
1226 call h5dcreate_f(this%active_group_id, trim(
field%name), precision_hdf, &
1227 filespace, dset_id, ierr, dcpl_id = dcpl_id)
1228 call h5sclose_f(filespace, ierr)
1229 call h5pclose_f(dcpl_id, ierr)
1235 dcount = [int(stride_ax_1, hsize_t), &
1236 int(stride_ax_2, hsize_t), &
1237 int(stride_ax_3, hsize_t), &
1238 int(counts, hsize_t)]
1239 doffset = [0_hsize_t, 0_hsize_t, 0_hsize_t, &
1240 int(offset, hsize_t) + append_offset]
1242 call h5dget_space_f(dset_id, filespace, ierr)
1244 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, doffset, dcount, &
1247 call h5screate_simple_f(dset_rank, dcount, memspace, ierr)
1252 if (this%precision == sp)
then
1257 call h5dwrite_f(dset_id, precision_hdf, write_buffer_sp, dcount, ierr, &
1258 file_space_id = filespace, mem_space_id = memspace, &
1260 deallocate(write_buffer_sp)
1261 else if (this%precision == dp)
then
1266 call h5dwrite_f(dset_id, precision_hdf, write_buffer_dp, dcount, ierr, &
1267 file_space_id = filespace, mem_space_id = memspace, &
1269 deallocate(write_buffer_dp)
1271 call neko_error(
"Unsupported precision")
1277 call h5pclose_f(xf_id, ierr)
1278 call h5sclose_f(memspace, ierr)
1279 call h5sclose_f(filespace, ierr)
1280 call h5dclose_f(dset_id, ierr)
1287 character(len=*),
intent(in) :: data_name
1288 type(vector_t),
intent(inout) :: vec
1289 character(len=*),
intent(in),
optional :: strategy
1290 character(len=1000) :: strategy_
1291 integer :: ierr, counts, offset, total_count, dset_rank
1292 integer(hid_t) :: precision_hdf
1293 integer(hid_t) :: xf_id, filespace, dset_id, memspace
1294 integer(hsize_t),
dimension(1) :: dcount, doffset
1295 integer(hsize_t),
dimension(1) :: tempddims, tempmaxddims
1297 logical :: dset_exists
1298 type(linear_dist_t) :: dist
1301 if (
present(strategy))
then
1302 if (trim(strategy) .eq.
"linear" .or. &
1303 trim(strategy) .eq.
"rank_0")
then
1304 strategy_ = strategy
1306 call neko_error(
"Unsupported strategy: " // trim(strategy))
1309 strategy_ =
"linear"
1318 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
1319 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
1320 precision_hdf = h5kind_to_type(rp, h5_real_kind)
1325 call h5lexists_f(this%active_group_id, trim(data_name), dset_exists, ierr)
1326 if (dset_exists)
then
1328 call h5dopen_f(this%active_group_id, trim(data_name), dset_id, ierr)
1330 call h5dget_space_f(dset_id, filespace, ierr)
1331 call h5sget_simple_extent_ndims_f(filespace, temprank, ierr)
1332 if (temprank .ne. 1)
then
1333 call neko_error(
"Dataset " // trim(data_name) // &
1334 " is not a rank 1 vector in file " // trim(
file_get_fname(this)))
1337 call h5sget_simple_extent_dims_f(filespace, tempddims, tempmaxddims, &
1339 call h5sclose_f(filespace, ierr)
1341 call neko_error(
"Dataset " // trim(data_name) // &
1342 " does not exist in current group " // trim(
file_get_fname(this)))
1348 total_count = int(tempddims(1))
1349 if (strategy_ .eq.
"linear")
then
1350 dist = linear_dist_t(total_count, pe_rank, pe_size, neko_comm)
1351 counts = dist%num_local()
1353 else if (strategy_ .eq.
"rank_0")
then
1354 if (pe_rank .eq. 0)
then
1355 counts = total_count
1361 call mpi_exscan(counts, offset, 1, mpi_integer, &
1362 mpi_sum, neko_comm, ierr)
1368 dcount = [int(counts, hsize_t)]
1369 doffset = [int(offset, hsize_t)]
1371 call h5dget_space_f(dset_id, filespace, ierr)
1373 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, doffset, dcount, &
1376 call h5screate_simple_f(dset_rank, dcount, memspace, ierr)
1381 call vec%init(counts, trim(data_name))
1382 call h5dread_f(dset_id, precision_hdf, vec%x, dcount, ierr, &
1383 file_space_id = filespace, mem_space_id = memspace, &
1389 call h5pclose_f(xf_id, ierr)
1390 call h5sclose_f(memspace, ierr)
1391 call h5sclose_f(filespace, ierr)
1392 call h5dclose_f(dset_id, ierr)
1399 character(len=*),
intent(in) :: data_name
1400 type(matrix_t),
intent(inout) :: mat
1401 character(len=*),
intent(in),
optional :: strategy
1402 character(len=1000) :: strategy_
1403 integer :: ierr, counts, offset, total_count, dset_rank
1404 integer(hid_t) :: precision_hdf
1405 integer(hid_t) :: xf_id, filespace, dset_id, memspace
1406 integer(hsize_t),
dimension(2) :: dcount, doffset
1407 integer(hsize_t),
dimension(2) :: tempddims, tempmaxddims
1409 logical :: dset_exists
1410 type(linear_dist_t) :: dist
1413 if (
present(strategy))
then
1414 if (trim(strategy) .eq.
"linear" .or. &
1415 trim(strategy) .eq.
"rank_0")
then
1416 strategy_ = strategy
1418 call neko_error(
"Unsupported strategy: " // trim(strategy))
1421 strategy_ =
"linear"
1430 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
1431 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
1432 precision_hdf = h5kind_to_type(rp, h5_real_kind)
1437 call h5lexists_f(this%active_group_id, trim(data_name), dset_exists, ierr)
1438 if (dset_exists)
then
1440 call h5dopen_f(this%active_group_id, trim(data_name), dset_id, ierr)
1442 call h5dget_space_f(dset_id, filespace, ierr)
1443 call h5sget_simple_extent_ndims_f(filespace, temprank, ierr)
1444 if (temprank .ne. 2)
then
1445 call neko_error(
"Dataset " // trim(data_name) // &
1446 " is not a rank 2 matrix in file " // trim(
file_get_fname(this)))
1449 call h5sget_simple_extent_dims_f(filespace, tempddims, tempmaxddims, &
1451 call h5sclose_f(filespace, ierr)
1453 call neko_error(
"Dataset " // trim(data_name) &
1454 //
" does not exist in current group " // &
1461 total_count = int(tempddims(2))
1462 if (strategy_ .eq.
"linear")
then
1463 dist = linear_dist_t(total_count, pe_rank, pe_size, neko_comm)
1464 counts = dist%num_local()
1466 else if (strategy_ .eq.
"rank_0")
then
1467 if (pe_rank .eq. 0)
then
1468 counts = total_count
1474 call mpi_scan(counts, offset, 1, mpi_integer, &
1475 mpi_sum, neko_comm, ierr)
1476 offset = offset - counts
1482 dcount = [int(tempddims(1), hsize_t), int(counts, hsize_t)]
1483 doffset = [0_hsize_t, int(offset, hsize_t)]
1485 call h5dget_space_f(dset_id, filespace, ierr)
1487 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, doffset, dcount, &
1490 call h5screate_simple_f(dset_rank, dcount, memspace, ierr)
1495 call mat%init(int(tempddims(1)), counts, trim(data_name))
1496 call h5dread_f(dset_id, precision_hdf, mat%x, dcount, ierr, &
1497 file_space_id = filespace, mem_space_id = memspace, &
1503 call h5pclose_f(xf_id, ierr)
1504 call h5sclose_f(memspace, ierr)
1505 call h5sclose_f(filespace, ierr)
1506 call h5dclose_f(dset_id, ierr)
1514 character(len=*),
intent(in) :: attr_name
1515 integer,
intent(in) :: attr
1517 integer(hid_t) :: filespace, attr_id
1518 integer(hsize_t),
dimension(1) :: dcount
1519 logical :: attr_exists
1524 dcount = [int(1, hsize_t)]
1525 call h5aexists_f(this%active_group_id, trim(attr_name), attr_exists, ierr)
1526 if (attr_exists)
then
1528 call h5aopen_f(this%active_group_id, trim(attr_name), attr_id, ierr)
1531 call h5screate_f(h5s_scalar_f, filespace, ierr)
1533 call h5acreate_f(this%active_group_id, trim(attr_name), &
1534 h5t_native_integer, &
1535 filespace, attr_id, ierr, h5p_default_f, h5p_default_f)
1536 call h5sclose_f(filespace, ierr)
1542 call h5awrite_f(attr_id, h5t_native_integer, attr, dcount, ierr)
1547 call h5aclose_f(attr_id, ierr)
1554 character(len=*),
intent(in) :: attr_name
1555 real(kind=rp),
intent(in) :: attr
1557 integer(hid_t) :: precision_hdf
1558 integer(hid_t) :: filespace, attr_id
1559 integer(hsize_t),
dimension(1) :: dcount
1560 logical :: attr_exists
1563 precision_hdf = h5kind_to_type(rp, h5_real_kind)
1568 dcount = [int(1, hsize_t)]
1569 call h5aexists_f(this%active_group_id, trim(attr_name), attr_exists, ierr)
1570 if (attr_exists)
then
1572 call h5aopen_f(this%active_group_id, trim(attr_name), attr_id, ierr)
1575 call h5screate_f(h5s_scalar_f, filespace, ierr)
1577 call h5acreate_f(this%active_group_id, trim(attr_name), precision_hdf, &
1578 filespace, attr_id, ierr, h5p_default_f, h5p_default_f)
1579 call h5sclose_f(filespace, ierr)
1585 call h5awrite_f(attr_id, precision_hdf, attr, dcount, ierr)
1590 call h5aclose_f(attr_id, ierr)
1597 character(len=*),
intent(in) :: attr_name
1598 integer,
intent(inout) :: attr
1599 logical,
intent(inout) :: attr_exists
1601 integer(hid_t) :: filespace, attr_id
1602 integer(hsize_t),
dimension(1) :: dcount
1607 dcount = [int(1, hsize_t)]
1608 call h5aexists_f(this%active_group_id, trim(attr_name), attr_exists, ierr)
1609 if (attr_exists)
then
1611 call h5aopen_f(this%active_group_id, trim(attr_name), attr_id, ierr)
1619 call h5aread_f(attr_id, h5t_native_integer, attr, dcount, ierr)
1624 call h5aclose_f(attr_id, ierr)
1631 character(len=*),
intent(in) :: attr_name
1632 real(kind=rp),
intent(inout) :: attr
1633 logical,
intent(inout) :: attr_exists
1635 integer(hid_t) :: precision_hdf
1636 integer(hid_t) :: filespace, attr_id
1637 integer(hsize_t),
dimension(1) :: dcount
1640 precision_hdf = h5kind_to_type(rp, h5_real_kind)
1645 dcount = [int(1, hsize_t)]
1646 call h5aexists_f(this%active_group_id, trim(attr_name), attr_exists, ierr)
1647 if (attr_exists)
then
1649 call h5aopen_f(this%active_group_id, trim(attr_name), attr_id, ierr)
1657 call h5aread_f(attr_id, precision_hdf, attr, dcount, ierr)
1662 call h5aclose_f(attr_id, ierr)
1671 character(len=1),
intent(in) :: mode
1672 call neko_error(
'Neko needs to be built with HDF5 support')
1678 call neko_error(
'Neko needs to be built with HDF5 support')
1684 character(len=*),
intent(in),
optional :: group_name_path
1685 call neko_error(
'Neko needs to be built with HDF5 support')
1691 class(*),
target,
intent(in) :: data
1692 real(kind=rp),
intent(in),
optional :: t
1693 call neko_error(
'Neko needs to be built with HDF5 support')
1699 class(*),
target,
intent(inout) :: data
1700 call neko_error(
'Neko needs to be built with HDF5 support')
1705 class(*),
intent(inout) :: data
1706 call neko_error(
'Neko needs to be built with HDF5 support')
1711 character(len=*),
intent(in) :: data_name
1712 class(*),
intent(inout) :: data
1713 character(len=*),
intent(in),
optional :: strategy
1714 call neko_error(
'Neko needs to be built with HDF5 support')
1719 character(len=*),
intent(in) :: data_name
1720 class(*),
intent(inout) :: data
1721 call neko_error(
'Neko needs to be built with HDF5 support')
1726 character(len=*),
intent(in) :: data_name
1727 class(*),
intent(inout) :: data
1728 logical,
intent(inout) :: exist
1729 call neko_error(
'Neko needs to be built with HDF5 support')
1734 type(vector_t),
intent(inout) :: vec
1735 call neko_error(
'Neko needs to be built with HDF5 support')
1740 type(matrix_t),
intent(inout) :: mat
1741 call neko_error(
'Neko needs to be built with HDF5 support')
1746 type(field_t),
intent(inout) :: fld
1747 call neko_error(
'Neko needs to be built with HDF5 support')
1752 character(len=*),
intent(in) :: data_name
1753 type(vector_t),
intent(inout) :: vec
1754 character(len=*),
intent(in),
optional :: strategy
1755 call neko_error(
'Neko needs to be built with HDF5 support')
1760 character(len=*),
intent(in) :: data_name
1761 type(matrix_t),
intent(inout) :: mat
1762 character(len=*),
intent(in),
optional :: strategy
1763 call neko_error(
'Neko needs to be built with HDF5 support')
1768 character(len=*),
intent(in) :: attr_name
1769 integer,
intent(in) :: attr
1770 call neko_error(
'Neko needs to be built with HDF5 support')
1775 character(len=*),
intent(in) :: attr_name
1776 real(kind=rp),
intent(in) :: attr
1777 call neko_error(
'Neko needs to be built with HDF5 support')
1782 character(len=*),
intent(in) :: attr_name
1783 integer,
intent(inout) :: attr
1784 logical,
intent(inout) :: attr_exists
1785 call neko_error(
'Neko needs to be built with HDF5 support')
1790 character(len=*),
intent(in) :: attr_name
1791 real(kind=rp),
intent(inout) :: attr
1792 logical,
intent(inout) :: attr_exists
1793 call neko_error(
'Neko needs to be built with HDF5 support')
integer, public pe_size
MPI size of communicator.
integer, public pe_rank
MPI rank.
type(mpi_comm), public neko_comm
MPI communicator.
Defines practical data distributions.
Defines a mapping of the degrees of freedom.
Contains the field_serties_t type.
Module for file I/O operations.
subroutine hdf5_file_open(this, mode)
Open a HDF5 file in a given mode.
character(len=1024) function file_get_fname(this)
Return the file name with the start counter.
subroutine hdf5_file_write_field(this, fld)
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_read_dataset(this, data_name, data, strategy)
subroutine hdf5_file_write_int_attribute(this, attr_name, attr)
subroutine hdf5_file_write_attribute(this, data_name, data)
subroutine hdf5_file_read_vector(this, data_name, vec, strategy)
subroutine hdf5_file_write_vector(this, vec)
subroutine hdf5_file_read_matrix(this, data_name, mat, strategy)
subroutine hdf5_file_read_int_attribute(this, attr_name, attr, attr_exists)
subroutine hdf5_file_write_rp_attribute(this, attr_name, attr)
subroutine hdf5_file_set_precision(this, precision)
Set the precision for the output (single or double)
subroutine hdf5_file_read_rp_attribute(this, attr_name, attr, attr_exists)
subroutine hdf5_file_read_attribute(this, data_name, data, exist)
subroutine hdf5_file_close(this)
Close the file.
subroutine hdf5_file_set_group(this, group_name_path)
Set the active group for HDF5 files.
subroutine hdf5_file_write_matrix(this, mat)
subroutine hdf5_file_write_dataset(this, data)
subroutine hdf5_file_set_overwrite(this, overwrite)
Set the overwrite flag for HDF5 files.
integer, parameter, public neko_log_debug
Debug log level.
type(log_t), public neko_log
Global log stream.
integer, parameter, public log_size
integer, parameter, public dp
integer, parameter, public sp
integer, parameter, public rp
Global precision used in computations.
subroutine, public filename_split(fname, path, name, suffix)
Extract file name components.
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
pure integer function, public filename_suffix_pos(fname)
Find position (in the string) of a filename's suffix.
Load-balanced linear distribution .
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.