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

Commit

Permalink
Added collect and distribute methods for coarrays
Browse files Browse the repository at this point in the history
  • Loading branch information
daminton committed Mar 28, 2023
1 parent 2b62ee1 commit 7a06aff
Show file tree
Hide file tree
Showing 3 changed files with 141 additions and 52 deletions.
93 changes: 54 additions & 39 deletions src/swiftest/swiftest_driver.f90
Original file line number Diff line number Diff line change
Expand Up @@ -81,17 +81,17 @@ program swiftest_driver
!$ write(param%display_unit,'(a)') ' ------------------'
!$ write(param%display_unit,'(a,i3,/)') ' Number of threads = ', nthreads
!$ if (param%log_output) write(*,'(a,i3)') ' OpenMP: Number of threads = ',nthreads
#ifdef COARRAY
if (this_image() == 1) then
write(param%display_unit,'(a)') ' Coarray parameters:'
write(param%display_unit,'(a)') ' -------------------'
write(param%display_unit,'(a,i3)') ' Number of images = ', num_images()
if (param%log_output) write(*,'(a,i3)') ' Coarray: Number of images = ',num_images()
#ifdef COARRAY
! Only execute file file I/O and reporting on image 1
if (this_image() == 1) then
write(param%display_unit,'(a)') ' Coarray parameters:'
write(param%display_unit,'(a)') ' -------------------'
write(param%display_unit,'(a,i3)') ' Number of images = ', num_images()
if (param%log_output) write(*,'(a,i3)') ' Coarray: Number of images = ',num_images()
#endif

call nbody_system%initialize(param)
call nbody_system%initialize(param)

associate (system_history => nbody_system%system_history)
! If this is a new run, compute energy initial conditions (if energy tracking is turned on) and write the initial conditions to file.
call nbody_system%display_run_information(param, integration_timer, phase="first")
if (param%lenergy) then
Expand All @@ -102,32 +102,40 @@ program swiftest_driver
end if
call nbody_system%conservation_report(param, lterminal=.true.)
end if
call system_history%take_snapshot(param,nbody_system)
call nbody_system%system_history%take_snapshot(param,nbody_system)
call nbody_system%dump(param)

do iloop = istart, nloops
!> Step the nbody_system forward in time
call integration_timer%start()
call nbody_system%step(param, nbody_system%t, dt)
call integration_timer%stop()

nbody_system%t = t0 + iloop * dt

!> Evaluate any discards or collisional outcomes
call nbody_system%discard(param)

!> If the loop counter is at the output cadence value, append the data file with a single frame
if (istep_out > 0) then
iout = iout + 1
if ((iout == istep) .or. (iloop == nloops)) then
iout = 0
idump = idump + 1
if (ltstretch) then
nout = nout + 1
istep = floor(istep_out * fstep_out**nout, kind=I4B)
end if

call system_history%take_snapshot(param,nbody_system)
#ifdef COARRAY
! Distribute test particles to the various images
call nbody_system%coarray_distribute()
end if ! this_image() == 1
#endif
do iloop = istart, nloops
!> Step the nbody_system forward in time
call integration_timer%start()
call nbody_system%step(param, nbody_system%t, dt)
call integration_timer%stop()

nbody_system%t = t0 + iloop * dt

!> Evaluate any discards or collisional outcomes
call nbody_system%discard(param)

!> If the loop counter is at the output cadence value, append the data file with a single frame
if (istep_out > 0) then
iout = iout + 1
if ((iout == istep) .or. (iloop == nloops)) then
iout = 0
idump = idump + 1
if (ltstretch) then
nout = nout + 1
istep = floor(istep_out * fstep_out**nout, kind=I4B)
end if
#ifdef COARRAY
if (this_image() == 1) then
call nbody_system%coarray_collect()
#endif
call nbody_system%system_history%take_snapshot(param,nbody_system)

if (idump == dump_cadence) then
idump = 0
Expand All @@ -138,20 +146,27 @@ program swiftest_driver
call nbody_system%display_run_information(param, integration_timer)
call integration_timer%reset()
if (param%lenergy) call nbody_system%conservation_report(param, lterminal=.true.)

end if
#ifdef COARRAY
call nbody_system%coarray_distribute()
end if
sync all
#endif
end if
end if

end do
! Dump any remaining history if it exists
end do
! Dump any remaining history if it exists
#ifdef COARRAY
if (this_image() == 1) then
#endif
call nbody_system%coarray_collect()
call nbody_system%dump(param)
call system_history%dump(param)
call nbody_system%system_history%dump(param)
call nbody_system%display_run_information(param, integration_timer, phase="last")

end associate
#ifdef COARRAY
end if ! this_image() == 1
#endif

end associate

call base_util_exit(SUCCESS)
Expand Down
29 changes: 22 additions & 7 deletions src/swiftest/swiftest_module.f90
Original file line number Diff line number Diff line change
Expand Up @@ -393,7 +393,11 @@ module swiftest
procedure :: read_in => swiftest_io_read_in_system !! Reads the initial conditions for an nbody system
procedure :: write_frame_system => swiftest_io_write_frame_system !! Write a frame of input data from file
procedure :: obl_pot => swiftest_obl_pot_system !! Compute the contribution to the total gravitational potential due solely to the oblateness of the central body
procedure :: dealloc => swiftest_util_dealloc_system !! Deallocates all allocatables and resets all values to defaults. Acts as a base for a finalizer
#ifdef COARRAY
procedure :: coarray_collect => swiftest_util_coarray_collect_system !! Collects test particles from distributed images into image #1
procedure :: coarray_distribute => swiftest_util_coarray_distribute_system !! Distributes test particles from image #1 out to all images
#endif
procedure :: dealloc => swiftest_util_dealloc_system !! Deallocates all allocatables and resets all values to defaults. Acts as a base for a finalizer
procedure :: get_energy_and_momentum => swiftest_util_get_energy_and_momentum_system !! Calculates the total nbody_system energy and momentum
procedure :: get_idvals => swiftest_util_get_idvalues_system !! Returns an array of all id values in use in the nbody_system
procedure :: rescale => swiftest_util_rescale_system !! Rescales the nbody_system into a new set of units
Expand Down Expand Up @@ -1184,6 +1188,18 @@ module subroutine swiftest_util_append_tp(self, source, lsource_mask)
logical, dimension(:), intent(in) :: lsource_mask !! Logical mask indicating which elements to append to
end subroutine swiftest_util_append_tp

#ifdef COARRAY
module subroutine swiftest_util_coarray_collect_system(self)
implicit none
class(swiftest_nbody_system), codimension[*], intent(inout) :: self
end subroutine swiftest_util_coarray_collect_system

module subroutine swiftest_util_coarray_distribute_system(self)
implicit none
class(swiftest_nbody_system), codimension[*], intent(inout) :: self
end subroutine swiftest_util_coarray_distribute_system
#endif

module subroutine swiftest_util_coord_b2h_pl(self, cb)
implicit none
class(swiftest_pl), intent(inout) :: self !! Swiftest massive body object
Expand Down Expand Up @@ -1659,19 +1675,18 @@ end subroutine swiftest_util_snapshot_save

module subroutine swiftest_util_snapshot_system(self, param, nbody_system, t, arg)
implicit none
class(swiftest_storage), intent(inout) :: self !! Swiftest storage object
class(swiftest_parameters), intent(inout) :: param !! Current run configuration parameters
class(swiftest_storage), intent(inout) :: self !! Swiftest storage object
class(swiftest_parameters), intent(inout) :: param !! Current run configuration parameters
#ifdef COARRAY
class(swiftest_nbody_system), intent(inout) :: nbody_system[*] !! Swiftest nbody system object to store
#else
class(swiftest_nbody_system), intent(inout) :: nbody_system !! Swiftest nbody system object to store
class(swiftest_nbody_system), intent(inout) :: nbody_system !! Swiftest nbody system object to store
#endif
real(DP), intent(in), optional :: t !! Time of snapshot if different from nbody_system time
character(*), intent(in), optional :: arg !! Optional argument (needed for extended storage type used in encounter snapshots)
real(DP), intent(in), optional :: t !! Time of snapshot if different from nbody_system time
character(*), intent(in), optional :: arg !! Optional argument (needed for extended storage type used in collision snapshots)
end subroutine swiftest_util_snapshot_system
end interface


interface swiftest_util_sort
pure module subroutine swiftest_util_sort_i4b(arr)
implicit none
Expand Down
71 changes: 65 additions & 6 deletions src/swiftest/swiftest_util.f90
Original file line number Diff line number Diff line change
Expand Up @@ -354,6 +354,65 @@ module subroutine swiftest_util_append_tp(self, source, lsource_mask)
return
end subroutine swiftest_util_append_tp

#ifdef COARRAY
module subroutine swiftest_util_coarray_collect_system(self)
!! author: David A. Minton
!!
!! Distributes test particles from image #1 out to all images.
implicit none
! Arguments
class(swiftest_nbody_system), codimension[*], intent(inout) :: self
! Internals
integer(I4B) :: i,j
integer(I4B), dimension(num_images()) :: ntp
class(swiftest_tp), allocatable :: tp_img

ntp(this_image()) = self%tp%nbody
sync all
if (this_image() == 1) then
do i = 2, num_images()
allocate(tp_img, source=self[i]%tp)
call self%tp%append(tp_img,lsource_mask=[(.true., j = 1, ntp(i))])
deallocate(tp_img)
end do
end if

return
end subroutine swiftest_util_coarray_collect_system


module subroutine swiftest_util_coarray_distribute_system(self)
!! author: David A. Minton
!!
!! Distributes test particles from image #1 out to all images.
implicit none
! Arguments
class(swiftest_nbody_system), codimension[*], intent(inout) :: self
! Internals
integer(I4B) :: i, istart, iend, ntot, num_per_image
class(swiftest_tp), allocatable :: tp_orig
logical, dimension(:), allocatable :: lspill_list

sync all
ntot = self[1]%tp%nbody
if (ntot == 0) return
allocate(lspill_list(ntot))
allocate(tp_orig, source=self[1]%tp)
if (allocated(self%tp)) deallocate(self%tp)
num_per_image = ntot / num_images()
istart = (this_image() - 1) * num_per_image + 1
if (this_image() == num_images()) then
iend = ntot
else
iend = this_image() * num_per_image
end if
lspill_list(:) = .false.
lspill_list(istart:iend) = .true.
call tp_orig%spill(self%tp,lspill_list(:), ldestructive=.false.)

return
end subroutine swiftest_util_coarray_distribute_system
#endif

module subroutine swiftest_util_coord_h2b_pl(self, cb)
!! author: David A. Minton
Expand Down Expand Up @@ -3065,15 +3124,15 @@ module subroutine swiftest_util_snapshot_system(self, param, nbody_system, t, ar
!! Takes a snapshot of the nbody_system for later file storage
implicit none
! Arguments
class(swiftest_storage), intent(inout) :: self !! Swiftest storage object
class(swiftest_parameters), intent(inout) :: param !! Current run configuration parameters
class(swiftest_storage), intent(inout) :: self !! Swiftest storage object
class(swiftest_parameters), intent(inout) :: param !! Current run configuration parameters
#ifdef COARRAY
class(swiftest_nbody_system), intent(inout) :: nbody_system[*] !! Swiftest nbody system object to store
class(swiftest_nbody_system), intent(inout) :: nbody_system[*] !! Swiftest nbody system object to store
#else
class(swiftest_nbody_system), intent(inout) :: nbody_system !! Swiftest nbody system object to store
class(swiftest_nbody_system), intent(inout) :: nbody_system !! Swiftest nbody system object to store
#endif
real(DP), intent(in), optional :: t !! Time of snapshot if different from nbody_system time
character(*), intent(in), optional :: arg !! Optional argument (needed for extended storage type used in collision snapshots)
real(DP), intent(in), optional :: t !! Time of snapshot if different from nbody_system time
character(*), intent(in), optional :: arg !! Optional argument (needed for extended storage type used in collision snapshots)
! Internals
class(swiftest_nbody_system), allocatable :: snapshot

Expand Down

0 comments on commit 7a06aff

Please sign in to comment.