diff --git a/sources/io.F90 b/sources/io.F90 index c91cbfc..f3002cb 100644 --- a/sources/io.F90 +++ b/sources/io.F90 @@ -3080,7 +3080,7 @@ module io #ifdef MPI use mesh , only : redistribute_blocks #endif /* MPI */ - use mpitools, only : nprocs, npmax, nproc + use mpitools, only : nprocs, nproc implicit none @@ -3088,7 +3088,7 @@ module io character(len=255) :: fname integer(hid_t) :: file_id, grp_id - integer :: nfiles, last_id, n + integer :: nfiles, last_id, n, i, nd, nr, nl, nu, il, iu logical :: flag real(kind=8) :: deinj @@ -3139,7 +3139,7 @@ module io ! larger or equal to the number of files, and when we have less processors than ! files ! - if (nproc < nfiles) then + if (nfiles <= nprocs .and. nproc < nfiles) then write(fname, "(a,'r',i6.6,'_',i5.5,'.h5')") trim(respath), nrest, nproc inquire(file=fname, exist=flag) @@ -3179,70 +3179,85 @@ module io call print_message(loc, "Could not close '" // trim(fname) // "'!") return end if + end if ! nproc < nfiles -#ifdef MPI -! if there are more files than processes, read the remaining files by -! the last process and redistribute blocks after each processed file, -! otherwise only redistribute blocks +! if there are more files than processes, divide the files equally between +! processes ! if (nprocs < nfiles) then - do n = nprocs, nfiles - 1 - call change_blocks_process(n, npmax) + nl = 0 + nd = nfiles / nprocs + nr = mod(nfiles, nprocs) + do n = 0, nprocs - 1 + if (n < nr) then + il = n * (nd + 1) + iu = il + nd + else + il = n * nd + nr + iu = il + nd - 1 + end if + do i = il, iu + call change_blocks_process(i, n) + end do + if (n == nproc) then + nl = il + nu = iu + end if + end do - if (nproc == npmax) then - write(fname, "(a,'r',i6.6,'_',i5.5,'.h5')") trim(respath), nrest, n - inquire(file=fname, exist=flag) - if (.not. flag) then - call print_message(loc, & - "Restart snapshot '" // trim(fname) // "' not found!") - status = 1 - return - end if + do n = nl, nu - call h5fopen_f(fname, H5F_ACC_RDONLY_F, file_id, status) - if (status /= 0) then - call print_message(loc, "Could not open '" // trim(fname) // "'!") - return - end if - - call restore_datablocks_h5(file_id, status) - if (status /= 0) & - call print_message(loc, & - "Could not restore datablocks from '" // trim(fname) // "'!") - - call h5gopen_f(file_id, 'attributes', grp_id, status) - if (status /= 0) then - call print_message(loc, "Could not open group 'attributes'!") - return - end if - call restore_attribute_h5(grp_id, 'einj', & - H5T_NATIVE_DOUBLE, 1, deinj, status) - if (status /= 0) & - call print_message(loc, "Could not get the injected energy!") - einj = einj + deinj - call h5gclose_f(grp_id, status) - if (status /= 0) & - call print_message(loc, "Could not close group 'attributes'!") - - call h5fclose_f(file_id, status) - if (status /= 0) then - call print_message(loc, "Could not close '" // trim(fname) // "'!") - return - end if + write(fname, "(a,'r',i6.6,'_',i5.5,'.h5')") trim(respath), nrest, n + inquire(file=fname, exist=flag) + if (.not. flag) then + call print_message(loc, & + "Restart snapshot '" // trim(fname) // "' not found!") + status = 1 + return end if - call redistribute_blocks(status) + call h5fopen_f(fname, H5F_ACC_RDONLY_F, file_id, status) + if (status /= 0) then + call print_message(loc, "Could not open '" // trim(fname) // "'!") + return + end if + + call restore_datablocks_h5(file_id, status) + if (status /= 0) & + call print_message(loc, & + "Could not restore datablocks from '" // trim(fname) // "'!") + + call h5gopen_f(file_id, 'attributes', grp_id, status) + if (status /= 0) then + call print_message(loc, "Could not open group 'attributes'!") + return + end if + call restore_attribute_h5(grp_id, 'einj', & + H5T_NATIVE_DOUBLE, 1, deinj, status) + if (status /= 0) & + call print_message(loc, "Could not get the injected energy!") + einj = einj + deinj + call h5gclose_f(grp_id, status) + if (status /= 0) & + call print_message(loc, "Could not close group 'attributes'!") + + call h5fclose_f(file_id, status) + if (status /= 0) then + call print_message(loc, "Could not close '" // trim(fname) // "'!") + return + end if end do - else - call redistribute_blocks(status) end if -#endif /* MPI */ if (allocated(block_array)) deallocate(block_array) +#ifdef MPI + call redistribute_blocks(status) +#endif /* MPI */ + !------------------------------------------------------------------------------- ! end subroutine read_restart_snapshot_h5