IO, DRIVER: Allow for restart with different number of ghost zones.

Now we can restart from restart snapshots with different number of ghost
zones. Also, remove the old restart snapshot format support from
read_datablocks_h5().

Signed-off-by: Grzegorz Kowal <grzegorz@amuncode.org>
This commit is contained in:
Grzegorz Kowal 2019-02-11 11:02:36 -02:00
parent 6d421e6efb
commit 601adc8037
2 changed files with 116 additions and 108 deletions

View File

@ -41,7 +41,7 @@ program amun
use boundaries , only : print_boundaries, boundary_variables use boundaries , only : print_boundaries, boundary_variables
use coordinates , only : initialize_coordinates, finalize_coordinates use coordinates , only : initialize_coordinates, finalize_coordinates
use coordinates , only : print_coordinates use coordinates , only : print_coordinates
use coordinates , only : nb => bcells use coordinates , only : nn => bcells
use domains , only : initialize_domains, finalize_domains use domains , only : initialize_domains, finalize_domains
use equations , only : initialize_equations, finalize_equations use equations , only : initialize_equations, finalize_equations
use equations , only : print_equations use equations , only : print_equations
@ -284,7 +284,6 @@ program amun
call read_snapshot_parameter("eqsys" , eqsys , iret) call read_snapshot_parameter("eqsys" , eqsys , iret)
call read_snapshot_parameter("eos" , eos , iret) call read_snapshot_parameter("eos" , eos , iret)
call read_snapshot_parameter("ncells" , ncells , iret) call read_snapshot_parameter("ncells" , ncells , iret)
call read_snapshot_parameter("nghosts", nghosts, iret)
call read_snapshot_parameter("maxlev" , toplev , iret) call read_snapshot_parameter("maxlev" , toplev , iret)
call read_snapshot_parameter("rdims" , bdims , iret) call read_snapshot_parameter("rdims" , bdims , iret)
call read_snapshot_parameter("xmin" , xmin , iret) call read_snapshot_parameter("xmin" , xmin , iret)
@ -300,7 +299,6 @@ program amun
call get_parameter("equation_system" , eqsys ) call get_parameter("equation_system" , eqsys )
call get_parameter("equation_of_state", eos ) call get_parameter("equation_of_state", eos )
call get_parameter("ncells" , ncells ) call get_parameter("ncells" , ncells )
call get_parameter("nghosts" , nghosts)
call get_parameter("maxlev" , toplev ) call get_parameter("maxlev" , toplev )
call get_parameter("xblocks" , bdims(1)) call get_parameter("xblocks" , bdims(1))
call get_parameter("yblocks" , bdims(2)) call get_parameter("yblocks" , bdims(2))
@ -317,6 +315,10 @@ program amun
#endif /* NDIMS == 3 */ #endif /* NDIMS == 3 */
end if end if
! get the number of ghost zones
!
call get_parameter("nghosts", nghosts)
! get the execution termination parameters ! get the execution termination parameters
! !
call get_parameter("nmax", nmax) call get_parameter("nmax", nmax)
@ -346,9 +348,9 @@ program amun
ymin, ymax, zmin, zmax, master, iret) ymin, ymax, zmin, zmax, master, iret)
if (iret > 0) go to 170 if (iret > 0) go to 170
#if NDIMS == 3 #if NDIMS == 3
call initialize_blocks((/ nv, nv, nb, nb, nb /), master, iret) call initialize_blocks((/ nv, nv, nn, nn, nn /), master, iret)
#else /* NDIMS == 3 */ #else /* NDIMS == 3 */
call initialize_blocks((/ nv, nv, nb, nb, 1 /), master, iret) call initialize_blocks((/ nv, nv, nn, nn, 1 /), master, iret)
#endif /* NDIMS == 3 */ #endif /* NDIMS == 3 */
if (iret > 0) go to 160 if (iret > 0) go to 160
call initialize_operators(master, iret) call initialize_operators(master, iret)
@ -417,6 +419,10 @@ program amun
! !
call build_leaf_list() call build_leaf_list()
! update boundaries
!
call boundary_variables(time, dtnext)
else else
! generate the initial mesh, refine that mesh to the desired level according to ! generate the initial mesh, refine that mesh to the desired level according to

View File

@ -1993,7 +1993,7 @@ module io
integer(hid_t) :: gid integer(hid_t) :: gid
integer :: ierr, l integer :: ierr, l
integer :: lndims, lmblocks, lnleafs, llast_id integer :: lndims, lmblocks, lnleafs, llast_id
integer :: lncells, lnghost, lnproc, lnseeds integer :: lncells, lnproc, lnseeds
! local pointers ! local pointers
! !
@ -2030,7 +2030,6 @@ module io
call read_attribute(gid, 'nleafs' , lnleafs ) call read_attribute(gid, 'nleafs' , lnleafs )
call read_attribute(gid, 'last_id', llast_id) call read_attribute(gid, 'last_id', llast_id)
call read_attribute(gid, 'ncells' , lncells ) call read_attribute(gid, 'ncells' , lncells )
call read_attribute(gid, 'nghosts', lnghost )
call read_attribute(gid, 'nseeds' , lnseeds ) call read_attribute(gid, 'nseeds' , lnseeds )
call read_attribute(gid, 'step' , step ) call read_attribute(gid, 'step' , step )
call read_attribute(gid, 'isnap' , isnap ) call read_attribute(gid, 'isnap' , isnap )
@ -2062,13 +2061,6 @@ module io
, "The block dimensions do not match!" , "The block dimensions do not match!"
end if end if
! check the number of ghost layers
!
if (lnghost /= nghosts) then
write(error_unit,"('[',a,']: ',a)") trim(loc) &
, "The number of ghost layers does not match!"
end if
! allocate space for seeds ! allocate space for seeds
! !
allocate(seeds(lnseeds)) allocate(seeds(lnseeds))
@ -2947,7 +2939,7 @@ module io
! !
use blocks , only : block_meta, block_data, list_data use blocks , only : block_meta, block_data, list_data
use blocks , only : append_datablock, link_blocks use blocks , only : append_datablock, link_blocks
use coordinates , only : nn => bcells use coordinates , only : nn => bcells, ng => nghosts
use equations , only : nv use equations , only : nv
use hdf5 , only : hid_t, hsize_t use hdf5 , only : hid_t, hsize_t
use hdf5 , only : h5gopen_f, h5gclose_f, h5lexists_f use hdf5 , only : h5gopen_f, h5gclose_f, h5lexists_f
@ -2971,7 +2963,7 @@ module io
character(len=16) :: bname character(len=16) :: bname
integer(hid_t) :: gid, bid integer(hid_t) :: gid, bid
integer(kind=4) :: l, i integer(kind=4) :: l, i
integer :: dblocks, ierr integer :: dblocks, ncells, nghosts, nc, nb, ne, status
! local arrays ! local arrays
! !
@ -2979,8 +2971,7 @@ module io
! local allocatable arrays ! local allocatable arrays
! !
integer(kind=4), dimension(:) , allocatable :: id real(kind=8), dimension(:,:,:,:,:), allocatable :: uu
real(kind=8) , dimension(:,:,:,:,:), allocatable :: uv, qv
! local parameters ! local parameters
! !
@ -2990,15 +2981,17 @@ module io
! !
! read the number of data blocks ! read the number of data blocks
! !
call h5gopen_f(fid, 'attributes', gid, ierr) call h5gopen_f(fid, 'attributes', gid, status)
if (ierr /= 0) then if (status /= 0) then
write(error_unit,"('[',a,']: ',a)") trim(loc) & write(error_unit,"('[',a,']: ',a)") trim(loc) &
, "Cannot open the attribute group!" , "Cannot open the attribute group!"
return return
end if end if
call read_attribute(gid, 'dblocks', dblocks) call read_attribute(gid, 'dblocks', dblocks)
call h5gclose_f(gid, ierr) call read_attribute(gid, 'ncells' , ncells )
if (ierr /= 0) then call read_attribute(gid, 'nghosts', nghosts)
call h5gclose_f(gid, status)
if (status /= 0) then
write(error_unit,"('[',a,']: ',a)") trim(loc) & write(error_unit,"('[',a,']: ',a)") trim(loc) &
, "Cannot close the attribute group!" , "Cannot close the attribute group!"
return return
@ -3010,73 +3003,42 @@ module io
! open the group 'datablocks' ! open the group 'datablocks'
! !
call h5gopen_f(fid, 'datablocks', gid, ierr) call h5gopen_f(fid, 'datablocks', gid, status)
if (ierr /= 0) then if (status /= 0) then
write(error_unit,"('[',a,']: ',a)") trim(loc) & write(error_unit,"('[',a,']: ',a)") trim(loc) &
, "Cannot open the data block group!" , "Cannot open the data block group!"
return else
! prepare the number of cells in the stored blocks and corresponding indices
!
nc = ncells + 2 * nghosts
if (nghosts >= ng) then
nb = 1 + (nghosts - ng)
ne = nc - (nghosts - ng)
else
nb = 1 + (ng - nghosts)
ne = nn - (ng - nghosts)
end if end if
! fill out dimensions dm(:) ! fill out dimensions
! !
dm(1) = dblocks dm(1) = dblocks
dm(2) = nv dm(2) = nv
dm(3) = nn dm(3) = nc
dm(4) = nn dm(4) = nc
#if NDIMS == 3 #if NDIMS == 3
dm(5) = nn dm(5) = nc
#endif /* NDIMS == 3 */ #endif /* NDIMS == 3 */
! check if the old data format is used, otherwise use the new one ! allocate arrays for input data
! !
call h5lexists_f(gid, 'meta', flag, ierr) allocate(uu(3,dm(2),dm(3),dm(4),dm(5)), stat = status)
! restart files are in the old format ! check if allocation was successful
! !
if (flag) then if (status == 0) then
! allocate arrays to read data ! iterate over data blocks
!
allocate(id(dm(1)))
allocate(uv(dm(1),dm(2),dm(3),dm(4),dm(5)))
allocate(qv(dm(1),dm(2),dm(3),dm(4),dm(5)))
! read array data from the HDF5 file
!
call read_array(gid, 'meta', dm(1:1), id(:) )
call read_array(gid, 'uvar', dm(1:5), uv(:,:,:,:,:))
call read_array(gid, 'qvar', dm(1:5), qv(:,:,:,:,:))
! iterate over data blocks, allocate them and fill out their fields
!
do l = 1, dm(1)
! allocate and append to the end of the list a new datablock
!
call append_datablock(pdata)
! associate the corresponding meta block with the current data block
!
call link_blocks(block_array(id(l))%ptr, pdata)
! fill out the array of conservative and primitive variables
!
pdata%u(:,:,:,:) = uv(l,:,:,:,:)
pdata%q(:,:,:,:) = qv(l,:,:,:,:)
end do ! l = 1, dm(1)
! deallocate allocatable arrays
!
if (allocated(id)) deallocate(id)
if (allocated(uv)) deallocate(uv)
if (allocated(qv)) deallocate(qv)
! restart files are in the new format
!
else ! flag
! iterate over data blocks, allocate them and fill out their fields
! !
do l = 1, dm(1) do l = 1, dm(1)
@ -3090,7 +3052,8 @@ module io
! open the group for the current block fields ! open the group for the current block fields
! !
call h5gopen_f(gid, bname, bid, ierr) call h5gopen_f(gid, bname, bid, status)
if (status == 0) then
! get the id of the linked meta block ! get the id of the linked meta block
! !
@ -3100,27 +3063,66 @@ module io
! !
call link_blocks(block_array(i)%ptr, pdata) call link_blocks(block_array(i)%ptr, pdata)
! fill out the array of primitive and conservative variables ! read the data
! !
call read_array(bid, 'pvar' , dm(2:5), pdata%q (:,:,:,:)) call read_array(bid, 'pvar' , dm(2:5), uu(1,:,:,:,:))
call read_array(bid, 'cvar0', dm(2:5), pdata%u0(:,:,:,:)) call read_array(bid, 'cvar0', dm(2:5), uu(2,:,:,:,:))
call read_array(bid, 'cvar1', dm(2:5), pdata%u1(:,:,:,:)) call read_array(bid, 'cvar1', dm(2:5), uu(3,:,:,:,:))
! fill out the block arrays taking into account the change of nghosts
!
if (nghosts >= ng) then
#if NDIMS == 3
pdata%q (:,:,:,:) = uu(1,:,nb:ne,nb:ne,nb:ne)
pdata%u0(:,:,:,:) = uu(2,:,nb:ne,nb:ne,nb:ne)
pdata%u1(:,:,:,:) = uu(3,:,nb:ne,nb:ne,nb:ne)
#else /* NDIMS == 3 */
pdata%q (:,:,:,:) = uu(1,:,nb:ne,nb:ne, : )
pdata%u0(:,:,:,:) = uu(2,:,nb:ne,nb:ne, : )
pdata%u1(:,:,:,:) = uu(3,:,nb:ne,nb:ne, : )
#endif /* NDIMS == 3 */
else
#if NDIMS == 3
pdata%q (:,nb:ne,nb:ne,nb:ne) = uu(1,:,:,:,:)
pdata%u0(:,nb:ne,nb:ne,nb:ne) = uu(2,:,:,:,:)
pdata%u1(:,nb:ne,nb:ne,nb:ne) = uu(3,:,:,:,:)
#else /* NDIMS == 3 */
pdata%q (:,nb:ne,nb:ne, : ) = uu(1,:,:,:,:)
pdata%u0(:,nb:ne,nb:ne, : ) = uu(2,:,:,:,:)
pdata%u1(:,nb:ne,nb:ne, : ) = uu(3,:,:,:,:)
#endif /* NDIMS == 3 */
end if
! close the current data block group ! close the current data block group
! !
call h5gclose_f(bid, ierr) call h5gclose_f(bid, status)
if (status /= 0) then
write(error_unit,"('[',a,']: ',a)") trim(loc) &
, "Cannot close group '", trim(bname), "'!"
end if
else
write(error_unit,"('[',a,']: ',a)") trim(loc) &
, "Cannot open group '", trim(bname), "'!"
end if
end do ! l = 1, dm(1) end do ! l = 1, dm(1)
end if ! flag else ! allocate
write(error_unit,"('[',a,']: ',a)") trim(loc) &
, "Cannot allocate temporary array!"
end if ! allocate
! deallocate allocatable arrays
!
if (allocated(uu)) deallocate(uu)
! close the data block group ! close the data block group
! !
call h5gclose_f(gid, ierr) call h5gclose_f(gid, status)
if (ierr /= 0) then if (status /= 0) then
write(error_unit,"('[',a,']: ',a)") trim(loc) & write(error_unit,"('[',a,']: ',a)") trim(loc) &
, "Cannot close the data block group!" , "Cannot close the data block group!"
return end if
end if end if
end if ! dblocks > 0 end if ! dblocks > 0