Skip to content

Commit

Permalink
Add option of task-local checkpointing.
Browse files Browse the repository at this point in the history
  • Loading branch information
p-costa committed Jan 29, 2024
1 parent fef5061 commit 024df40
Show file tree
Hide file tree
Showing 2 changed files with 75 additions and 1 deletion.
74 changes: 74 additions & 0 deletions src/load.f90
Original file line number Diff line number Diff line change
Expand Up @@ -467,6 +467,80 @@ subroutine load_one(io,filename,comm,ng,nh,lo,hi,p,time,istep)
end select
end subroutine load_one
!
subroutine load_all_local(io,filename,n,nh,u,v,w,p,time,istep)
!
! reads/writes a restart file
!
implicit none
character(len=1), intent(in) :: io
character(len=*), intent(in) :: filename
integer , intent(in), dimension(3) :: n,nh
real(rp), intent(inout), dimension(1-nh(1):,1-nh(2):,1-nh(3):) :: u,v,w,p
real(rp), intent(inout) :: time
integer , intent(inout) :: istep
real(rp), dimension(2) :: fldinfo
integer :: iunit
!
select case(io)
case('r')
open(newunit=iunit,file=filename,action='read' ,access='stream',form='unformatted',status='old' )
read(iunit) u(1:n(1),1:n(2),1:n(3)), &
v(1:n(1),1:n(2),1:n(3)), &
w(1:n(1),1:n(2),1:n(3)), &
p(1:n(1),1:n(2),1:n(3)), &
fldinfo(1:2)
close(iunit)
time = fldinfo(1)
istep = nint(fldinfo(2))
case('w')
!
! write
!
fldinfo = [time,1._rp*istep]
open(newunit=iunit,file=filename,action='write',access='stream',form='unformatted',status='replace')
write(iunit) u(1:n(1),1:n(2),1:n(3)), &
v(1:n(1),1:n(2),1:n(3)), &
w(1:n(1),1:n(2),1:n(3)), &
p(1:n(1),1:n(2),1:n(3)), &
fldinfo(1:2)
close(iunit)
end select
end subroutine load_all_local
!
subroutine load_one_local(io,filename,n,nh,p,time,istep)
!
! reads/writes a restart file
!
implicit none
character(len=1), intent(in) :: io
character(len=*), intent(in) :: filename
integer , intent(in), dimension(3) :: n,nh
real(rp), intent(inout), dimension(1-nh(1):,1-nh(2):,1-nh(3):) :: p
real(rp), intent(inout), optional :: time
integer , intent(inout), optional :: istep
real(rp), dimension(2) :: fldinfo
integer :: iunit
!
select case(io)
case('r')
open(newunit=iunit,file=filename,action='read' ,access='stream',form='unformatted',status='old' )
read(iunit) p(1:n(1),1:n(2),1:n(3))
if(present(time) .and. present(istep)) read(iunit) fldinfo(1:2)
close(iunit)
time = fldinfo(1)
istep = nint(fldinfo(2))
case('w')
!
! write
!
fldinfo = [time,1._rp*istep]
open(newunit=iunit,file=filename,action='write',access='stream',form='unformatted',status='replace')
write(iunit) p(1:n(1),1:n(2),1:n(3))
if(present(time) .and. present(istep)) write(iunit) fldinfo(1:2)
close(iunit)
end select
end subroutine load_one_local
!
#if defined(_USE_HDF5)
subroutine io_field_hdf5(io,filename,varname,ng,nh,lo,hi,var,meta,x_g,y_g,z_g)
!
Expand Down
2 changes: 1 addition & 1 deletion src/main.f90
Original file line number Diff line number Diff line change
Expand Up @@ -213,7 +213,7 @@ program cans
if(myid == 0) print*, ''
call initgrid(gtype,ng(3),gr,l(3),dzc_g,dzf_g,zc_g,zf_g)
if(myid == 0) then
open(99,file=trim(datadir)//'grid.bin',access='stream')
open(99,file=trim(datadir)//'grid.bin',action='write',format='unformatted',access='stream',status='replace')
write(99) dzc_g(1:ng(3)),dzf_g(1:ng(3)),zc_g(1:ng(3)),zf_g(1:ng(3))
close(99)
open(99,file=trim(datadir)//'grid.out')
Expand Down

0 comments on commit 024df40

Please sign in to comment.