Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add option of task-local checkpointing. #110

Merged
merged 2 commits into from
Mar 5, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
82 changes: 78 additions & 4 deletions src/load.f90
Original file line number Diff line number Diff line change
Expand Up @@ -467,8 +467,83 @@ 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)
use hdf5
!
! collective single field data I/O using HDF5
!
Expand Down Expand Up @@ -501,7 +576,6 @@ subroutine io_field_hdf5(io,filename,varname,ng,nh,lo,hi,var,meta,x_g,y_g,z_g)
integer(HSIZE_T) , dimension(3) :: data_count
integer(HSSIZE_T), dimension(3) :: data_offset
integer(HSSIZE_T), dimension(3) :: halo_offset
type(MPI_INFO) :: info = MPI_INFO_NULL
!
n(:) = hi(:)-lo(:)+1
sizes(:) = ng(:)
Expand All @@ -517,7 +591,7 @@ subroutine io_field_hdf5(io,filename,varname,ng,nh,lo,hi,var,meta,x_g,y_g,z_g)
select case(io)
case('r')
call h5pcreate_f(H5P_FILE_ACCESS_F,plist_id,ierr)
call h5pset_fapl_mpio_f(plist_id,MPI_COMM_WORLD%MPI_VAL,info%MPI_VAL,ierr)
call h5pset_fapl_mpio_f(plist_id,MPI_COMM_WORLD,MPI_INFO_NULL,ierr)
call h5fopen_f(filename,H5F_ACC_RDONLY_F,file_id,ierr,access_prp=plist_id)
call h5pclose_f(plist_id,ierr)
!
Expand All @@ -544,11 +618,11 @@ subroutine io_field_hdf5(io,filename,varname,ng,nh,lo,hi,var,meta,x_g,y_g,z_g)
call h5dclose_f(dset,ierr)
call h5fclose_f(file_id,ierr)
end if
call MPI_Bcast(meta,2,MPI_REAL_RP,0,MPI_COMM_WORLD)
call MPI_Bcast(meta,2,MPI_REAL_RP,0,MPI_COMM_WORLD,ierr)
case('w')
call h5screate_simple_f(ndims,dims,filespace,ierr)
call h5pcreate_f(H5P_FILE_ACCESS_F,plist_id,ierr)
call h5pset_fapl_mpio_f(plist_id,MPI_COMM_WORLD%MPI_VAL,info%MPI_VAL,ierr)
call h5pset_fapl_mpio_f(plist_id,MPI_COMM_WORLD,MPI_INFO_NULL,ierr)
call h5fcreate_f(filename,H5F_ACC_TRUNC_F,file_id,ierr,access_prp=plist_id)
call h5pclose_f(plist_id,ierr)
!
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',form='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
Loading