Skip to content
This repository was archived by the owner on Aug 28, 2024. It is now read-only.

Commit

Permalink
Fixed restarts so that data at a specific time value are read in
Browse files Browse the repository at this point in the history
  • Loading branch information
daminton committed Jan 6, 2023
1 parent 6a3e4b5 commit 197347d
Show file tree
Hide file tree
Showing 5 changed files with 52 additions and 22 deletions.
36 changes: 33 additions & 3 deletions src/netcdf_io/netcdf_io_implementations.f90
Original file line number Diff line number Diff line change
Expand Up @@ -15,7 +15,6 @@ module subroutine netcdf_io_check(status, call_identifier)
!! author: Carlisle A. Wishard, Dana Singh, and David A. Minton
!!
!! Checks the status of all NetCDF operations to catch errors
use netcdf
implicit none
! Arguments
integer, intent (in) :: status !! The status code returned by a NetCDF function
Expand All @@ -35,7 +34,6 @@ module subroutine netcdf_io_close(self)
!! author: Carlisle A. Wishard, Dana Singh, and David A. Minton
!!
!! Closes a NetCDF file
use netcdf
implicit none
! Arguments
class(netcdf_parameters), intent(inout) :: self !! Parameters used to identify a particular NetCDF dataset
Expand All @@ -49,12 +47,44 @@ module subroutine netcdf_io_close(self)
end subroutine netcdf_io_close



module subroutine netcdf_io_find_tslot(self, t)
!! author: David A. Minton
!!
!! Given an open NetCDF file and a value of time t, finds the index of the time value (aka the time slot) to place a new set of data.
!! The returned value of tslot will correspond to the first index value where the value of t is greater than or equal to the saved time value.
implicit none
! Arguments
class(netcdf_parameters), intent(inout) :: self !! Parameters used to identify a particular NetCDF dataset
real(DP), intent(in) :: t !! The value of time to search for
! Internals
real(DP), dimension(:), allocatable :: tvals

self%tslot = 0

if (.not.self%lfile_is_open) return

call netcdf_io_check( nf90_inquire_dimension(self%id, self%time_dimid, self%time_dimname, len=self%max_tslot), "netcdf_io_find_tslot nf90_inquire_dimension max_tslot" )
if (self%max_tslot > 0) then
allocate(tvals(self%max_tslot))
call netcdf_io_check( nf90_get_var(self%id, self%time_varid, tvals(:), start=[1]), "netcdf_io_find_tslot get_var" )
else
allocate(tvals(1))
tvals(1) = -huge(1.0_DP)
end if

self%tslot = findloc(tvals, t, dim=1)
if (self%tslot == 0) self%tslot = self%max_tslot + 1
self%max_tslot = max(self%max_tslot, self%tslot)

return
end subroutine netcdf_io_find_tslot

module subroutine netcdf_io_sync(self)
!! author: David A. Minton
!!
!! Syncrhonize the disk and memory buffer of the NetCDF file (e.g. commit the frame files stored in memory to disk)
!!
use netcdf
implicit none
! Arguments
class(netcdf_parameters), intent(inout) :: self !! Parameters used to identify a particular NetCDF dataset
Expand Down
11 changes: 9 additions & 2 deletions src/netcdf_io/netcdf_io_module.f90
Original file line number Diff line number Diff line change
Expand Up @@ -138,8 +138,9 @@ module netcdf_io
integer(I4B) :: discard_body_id_varid !! ID for the id of the other body involved in the discard
logical :: lpseudo_vel_exists = .false. !! Logical flag to indicate whether or not the pseudovelocity vectors were present in an old file.
contains
procedure :: close => netcdf_io_close !! Closes an open NetCDF file
procedure :: sync => netcdf_io_sync !! Syncrhonize the disk and memory buffer of the NetCDF file (e.g. commit the frame files stored in memory to disk)
procedure :: close => netcdf_io_close !! Closes an open NetCDF file
procedure :: find_tslot => netcdf_io_find_tslot !! Finds the time dimension index for a given value of t
procedure :: sync => netcdf_io_sync !! Syncrhonize the disk and memory buffer of the NetCDF file (e.g. commit the frame files stored in memory to disk)
end type netcdf_parameters

interface
Expand All @@ -153,6 +154,12 @@ module subroutine netcdf_io_close(self)
implicit none
class(netcdf_parameters), intent(inout) :: self !! Parameters used to identify a particular NetCDF dataset
end subroutine netcdf_io_close

module subroutine netcdf_io_find_tslot(self, t)
implicit none
class(netcdf_parameters), intent(inout) :: self !! Parameters used to identify a particular NetCDF dataset
real(DP), intent(in) :: t !! The value of time to search for
end subroutine netcdf_io_find_tslot

module subroutine netcdf_io_sync(self)
implicit none
Expand Down
6 changes: 3 additions & 3 deletions src/swiftest/swiftest_driver.f90
Original file line number Diff line number Diff line change
Expand Up @@ -63,6 +63,8 @@ program swiftest_driver
nloops = ceiling((tstop - t0) / dt, kind=I8B)
istart = ceiling((tstart - t0) / dt + 1.0_DP, kind=I8B)
iloop = istart - 1
iout = 0
idump = 0

! Set up nbody_system storage for intermittent file dumps
if (dump_cadence == 0) dump_cadence = ceiling(nloops / (1.0_DP * istep_out), kind=I8B)
Expand Down Expand Up @@ -104,9 +106,7 @@ program swiftest_driver
call nbody_system%compact_output(param,integration_timer)
end if

iout = 0
idump = 0
nbody_system%t = tstart

do iloop = istart, nloops
!> Step the nbody_system forward in time
call integration_timer%start()
Expand Down
19 changes: 5 additions & 14 deletions src/swiftest/swiftest_io.f90
Original file line number Diff line number Diff line change
Expand Up @@ -278,27 +278,16 @@ module subroutine swiftest_io_dump_storage(self, param)
class(swiftest_parameters), intent(inout) :: param !! Current run configuration parameters
! Internals
integer(I4B) :: i
real(DP), dimension(:), allocatable :: tvals

if (self%iframe == 0) return
call self%make_index_map()
associate(nc => param%system_history%nc)
call nc%open(param)
! Get the current time values in the file if this is an old file
if (nc%max_tslot > 0) then
allocate(tvals(nc%max_tslot))
call netcdf_io_check( nf90_get_var(nc%id, nc%time_varid, tvals(:), start=[1]), "swiftest_io_dump_storage time_varid" )
else
allocate(tvals(1))
tvals(1) = -huge(1.0_DP)
end if

do i = 1, self%iframe
if (allocated(self%frame(i)%item)) then
select type(nbody_system => self%frame(i)%item)
class is (swiftest_nbody_system)
nc%tslot = findloc(tvals, nbody_system%t, dim=1)
if (nc%tslot == 0) nc%tslot = nc%max_tslot + 1
nc%max_tslot = max(nc%max_tslot, nc%tslot)
call nbody_system%write_frame(param)
end select
deallocate(self%frame(i)%item)
Expand Down Expand Up @@ -907,6 +896,7 @@ module function swiftest_io_netcdf_read_frame_system(self, nc, param) result(ier


call nc%open(param, readonly=.true.)
call nc%find_tslot(self%t)
call self%read_hdr(nc, param)
associate(cb => self%cb, pl => self%pl, tp => self%tp, npl => self%pl%nbody, ntp => self%tp%nbody, tslot => nc%tslot)

Expand Down Expand Up @@ -1585,10 +1575,11 @@ module subroutine swiftest_io_netcdf_write_frame_system(self, nc, param)
!! Write a frame (header plus records for each massive body and active test particle) to a output binary file
implicit none
! Arguments
class(swiftest_nbody_system), intent(inout) :: self !! Swiftest nbody_system object
class(swiftest_nbody_system), intent(inout) :: self !! Swiftest nbody_system object
class(swiftest_netcdf_parameters), intent(inout) :: nc !! Parameters used to for writing a NetCDF dataset to file
class(swiftest_parameters), intent(inout) :: param !! Current run configuration parameters
class(swiftest_parameters), intent(inout) :: param !! Current run configuration parameters

call nc%find_tslot(self%t)
call self%write_hdr(nc, param)
call self%cb%write_frame(nc, param)
call self%pl%write_frame(nc, param)
Expand Down
2 changes: 2 additions & 0 deletions src/swiftest/swiftest_util.f90
Original file line number Diff line number Diff line change
Expand Up @@ -2645,6 +2645,8 @@ module subroutine swiftest_util_setup_construct_system(nbody_system, param)
end select
call nbody_system%collider%setup(nbody_system)

nbody_system%t = param%tstart


return
end subroutine swiftest_util_setup_construct_system
Expand Down

0 comments on commit 197347d

Please sign in to comment.