IO: Resume jobs from HDF5 restart snapshots for different nghosts.
Signed-off-by: Grzegorz Kowal <grzegorz@amuncode.org>
This commit is contained in:
parent
1f1ae8084f
commit
a03375122b
130
sources/io.F90
130
sources/io.F90
@ -3613,7 +3613,7 @@ module io
|
||||
subroutine store_attributes_h5(loc_id, problem, restart, status)
|
||||
|
||||
use blocks , only : get_mblocks, get_dblocks, get_nleafs
|
||||
use blocks , only : get_last_id, nchildren
|
||||
use blocks , only : get_last_id, nregs, nchildren
|
||||
use coordinates, only : minlev, maxlev
|
||||
use coordinates, only : bcells, ncells, nghosts
|
||||
use coordinates, only : bdims => domain_base_dims
|
||||
@ -3730,6 +3730,8 @@ module io
|
||||
H5T_NATIVE_INTEGER, 1, nchildren, status)
|
||||
call store_attribute_h5(grp_id, 'mblocks', &
|
||||
H5T_NATIVE_INTEGER, 1, get_mblocks(), status)
|
||||
call store_attribute_h5(grp_id, 'nregisters', &
|
||||
H5T_NATIVE_INTEGER, 1, nregs, status)
|
||||
call store_attribute_h5(grp_id, 'last_id', &
|
||||
H5T_NATIVE_INTEGER, 1, get_last_id(), status)
|
||||
call store_attribute_h5(grp_id, 'maximum_rejections', &
|
||||
@ -4253,7 +4255,8 @@ module io
|
||||
subroutine restore_datablocks_h5(loc_id, status)
|
||||
|
||||
use blocks , only : block_meta, block_data
|
||||
use blocks , only : append_datablock, link_blocks
|
||||
use blocks , only : append_datablock, link_blocks, nregs
|
||||
use coordinates, only : bcells, nghosts
|
||||
use helpers , only : print_message
|
||||
|
||||
implicit none
|
||||
@ -4265,10 +4268,11 @@ module io
|
||||
|
||||
integer(hid_t) :: grp_id, blk_id
|
||||
character(len=64) :: blk_name
|
||||
integer(kind=4) :: dblocks, l, id
|
||||
integer(kind=4) :: dblocks, l, id, nr, nv, nm, ng, nl, nu
|
||||
|
||||
integer(hsize_t), dimension(4) :: pdims = 1
|
||||
integer(hsize_t), dimension(5) :: cdims = 1
|
||||
integer(hsize_t), dimension(5) :: dims = 1
|
||||
|
||||
real(kind=8), dimension(:,:,:,:,:), allocatable :: array
|
||||
|
||||
character(len=*), parameter :: loc = 'IO::restore_datablocks_h5()'
|
||||
|
||||
@ -4282,6 +4286,14 @@ module io
|
||||
|
||||
call restore_attribute_h5(grp_id, 'dblocks', &
|
||||
H5T_NATIVE_INTEGER, 1, dblocks, status)
|
||||
call restore_attribute_h5(grp_id, 'nregisters', &
|
||||
H5T_NATIVE_INTEGER, 1, nr, status)
|
||||
call restore_attribute_h5(grp_id, 'nvars', &
|
||||
H5T_NATIVE_INTEGER, 1, nv, status)
|
||||
call restore_attribute_h5(grp_id, 'bcells', &
|
||||
H5T_NATIVE_INTEGER, 1, nm, status)
|
||||
call restore_attribute_h5(grp_id, 'nghosts', &
|
||||
H5T_NATIVE_INTEGER, 1, ng, status)
|
||||
|
||||
call h5gclose_f(grp_id, status)
|
||||
if (status /= 0) then
|
||||
@ -4297,39 +4309,89 @@ module io
|
||||
return
|
||||
end if
|
||||
|
||||
do l = 1, dblocks
|
||||
write(blk_name, "('datablock_', i0)") l
|
||||
#if NDIMS == 3
|
||||
allocate(array(nv,nm,nm,nm,nr), stat=status)
|
||||
#else /* NDIMS == 3 */
|
||||
allocate(array(nv,nm,nm, 1,nr), stat=status)
|
||||
#endif /* NDIMS == 3 */
|
||||
if (status == 0) then
|
||||
dims = shape(array)
|
||||
|
||||
call append_datablock(pdata, status)
|
||||
if (status /= 0) then
|
||||
call print_message(loc, "Could not append new datablock!")
|
||||
go to 1000
|
||||
end if
|
||||
nr = min(nr, nregs)
|
||||
|
||||
call h5gopen_f(grp_id, blk_name, blk_id, status)
|
||||
if (status == 0) then
|
||||
call restore_attribute_h5(blk_id, 'meta', &
|
||||
H5T_NATIVE_INTEGER, 1, id, status)
|
||||
|
||||
call link_blocks(block_array(id)%ptr, pdata)
|
||||
|
||||
pdims = shape(pdata%q)
|
||||
cdims = shape(pdata%uu)
|
||||
|
||||
call read_dataset_h5(blk_id, 'primitive_variables', &
|
||||
H5T_NATIVE_DOUBLE, pdims, pdata%q , status)
|
||||
call read_dataset_h5(blk_id, 'conservative_variables', &
|
||||
H5T_NATIVE_DOUBLE, cdims, pdata%uu, status)
|
||||
|
||||
call h5gclose_f(blk_id, status)
|
||||
if (status /= 0) &
|
||||
call print_message(loc, &
|
||||
"Could not close group '" // trim(blk_name) // "'!")
|
||||
if (ng >= nghosts) then
|
||||
nl = 1 + (ng - nghosts)
|
||||
nu = nm - (ng - nghosts)
|
||||
else
|
||||
call print_message(loc, &
|
||||
"Could not open group '" // trim(blk_name) // "'!")
|
||||
nl = 1 + (nghosts - ng)
|
||||
nu = bcells - (nghosts - ng)
|
||||
end if
|
||||
end do
|
||||
|
||||
do l = 1, dblocks
|
||||
write(blk_name, "('datablock_', i0)") l
|
||||
|
||||
call append_datablock(pdata, status)
|
||||
if (status /= 0) then
|
||||
call print_message(loc, "Could not append new datablock!")
|
||||
go to 1000
|
||||
end if
|
||||
|
||||
call h5gopen_f(grp_id, blk_name, blk_id, status)
|
||||
if (status == 0) then
|
||||
call restore_attribute_h5(blk_id, 'meta', &
|
||||
H5T_NATIVE_INTEGER, 1, id, status)
|
||||
|
||||
call link_blocks(block_array(id)%ptr, pdata)
|
||||
|
||||
call read_dataset_h5(blk_id, 'primitive_variables', &
|
||||
H5T_NATIVE_DOUBLE, dims(1:4), array(:,:,:,:,1), status)
|
||||
if (ng >= nghosts) then
|
||||
#if NDIMS == 3
|
||||
pdata%q(:,:,:,:) = array(:,nl:nu,nl:nu,nl:nu,1)
|
||||
#else /* NDIMS == 3 */
|
||||
pdata%q(:,:,:,:) = array(:,nl:nu,nl:nu, : ,1)
|
||||
#endif /* NDIMS == 3 */
|
||||
else
|
||||
#if NDIMS == 3
|
||||
pdata%q(:,nl:nu,nl:nu,nl:nu) = array(:,:,:,:,1)
|
||||
#else /* NDIMS == 3 */
|
||||
pdata%q(:,nl:nu,nl:nu, : ) = array(:,:,:,:,1)
|
||||
#endif /* NDIMS == 3 */
|
||||
end if
|
||||
call read_dataset_h5(blk_id, 'conservative_variables', &
|
||||
H5T_NATIVE_DOUBLE, dims(1:5), array(:,:,:,:,:), status)
|
||||
if (ng >= nghosts) then
|
||||
#if NDIMS == 3
|
||||
pdata%uu(:,:,:,:,1:nr) = array(:,nl:nu,nl:nu,nl:nu,1:nr)
|
||||
#else /* NDIMS == 3 */
|
||||
pdata%uu(:,:,:,:,1:nr) = array(:,nl:nu,nl:nu, : ,1:nr)
|
||||
#endif /* NDIMS == 3 */
|
||||
else
|
||||
#if NDIMS == 3
|
||||
pdata%uu(:,nl:nu,nl:nu,nl:nu,1:nr) = array(:,:,:,:,1:nr)
|
||||
#else /* NDIMS == 3 */
|
||||
pdata%uu(:,nl:nu,nl:nu, : ,1:nr) = array(:,:,:,:,1:nr)
|
||||
#endif /* NDIMS == 3 */
|
||||
end if
|
||||
|
||||
call h5gclose_f(blk_id, status)
|
||||
if (status /= 0) &
|
||||
call print_message(loc, &
|
||||
"Could not close group '" // trim(blk_name) // "'!")
|
||||
else
|
||||
call print_message(loc, &
|
||||
"Could not open group '" // trim(blk_name) // "'!")
|
||||
end if
|
||||
end do
|
||||
|
||||
deallocate(array, stat=status)
|
||||
if (status /= 0) &
|
||||
call print_message(loc, &
|
||||
"Could not deallocate memory for the temporary arrays!")
|
||||
else
|
||||
call print_message(loc, &
|
||||
"Could not allocate memory for the temporary arrays!")
|
||||
end if
|
||||
|
||||
1000 continue
|
||||
|
||||
|
Loading…
x
Reference in New Issue
Block a user