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

Commit

Permalink
Removed spurious move_alloc that was accidentally copying over the ol…
Browse files Browse the repository at this point in the history
…d encounter list (with incorrect indices) over to the new system. Also added an additional symba_util_copy_encounter to ensure all quantities in an encounter list get copied over during rearray
  • Loading branch information
daminton committed Aug 17, 2021
1 parent d44b452 commit b6d35fb
Show file tree
Hide file tree
Showing 3 changed files with 38 additions and 9 deletions.
14 changes: 7 additions & 7 deletions src/io/io.f90
Original file line number Diff line number Diff line change
Expand Up @@ -72,13 +72,13 @@ module subroutine io_conservation_report(self, param, lterminal)
Etotal_error = (Eorbit_now - Ecollisions - Eorbit_orig - Euntracked) / abs(Eorbit_orig)
Merror = (Mtot_now - Mtot_orig) / Mtot_orig
write(*, egytermfmt) Lerror, Ecoll_error, Etotal_error, Merror
if (Ecoll_error > 0.0_DP) then
write(*,*) 'Something has gone wrong! Collisional energy should not be positive!'
write(*,*) 'dke_orbit: ',(ke_orbit_now - ke_orbit_last) / abs(Eorbit_orig)
write(*,*) 'dke_spin : ',(ke_spin_now - ke_spin_last) / abs(Eorbit_orig)
write(*,*) 'dpe : ',(pe_now - pe_last) / abs(Eorbit_orig)
write(*,*)
end if
! if (Ecoll_error > 0.0_DP) then
! write(*,*) 'Something has gone wrong! Collisional energy should not be positive!'
! write(*,*) 'dke_orbit: ',(ke_orbit_now - ke_orbit_last) / abs(Eorbit_orig)
! write(*,*) 'dke_spin : ',(ke_spin_now - ke_spin_last) / abs(Eorbit_orig)
! write(*,*) 'dpe : ',(pe_now - pe_last) / abs(Eorbit_orig)
! write(*,*)
! end if
end if
ke_orbit_last = ke_orbit_now
ke_spin_last = ke_spin_now
Expand Down
8 changes: 8 additions & 0 deletions src/modules/symba_classes.f90
Original file line number Diff line number Diff line change
Expand Up @@ -140,6 +140,7 @@ module symba_classes
procedure :: encounter_check => symba_encounter_check !! Checks if massive bodies are going through close encounters with each other
procedure :: kick => symba_kick_encounter !! Kick barycentric velocities of active test particles within SyMBA recursion
procedure :: setup => symba_setup_encounter !! A constructor that sets the number of encounters and allocates and initializes all arrays
procedure :: copy => symba_util_copy_encounter !! Copies elements from the source encounter list into self.
procedure :: spill => symba_util_spill_encounter !! "Spills" bodies from one object to another depending on the results of a mask (uses the PACK intrinsic)
procedure :: append => symba_util_append_encounter !! Appends elements from one structure to another
end type symba_encounter
Expand Down Expand Up @@ -551,6 +552,13 @@ module subroutine symba_util_append_tp(self, source, lsource_mask)
class(swiftest_body), intent(in) :: source !! Source object to append
logical, dimension(:), intent(in) :: lsource_mask !! Logical mask indicating which elements to append to
end subroutine symba_util_append_tp

module subroutine symba_util_copy_encounter(self, source)
use swiftest_classes, only : swiftest_encounter
implicit none
class(symba_encounter), intent(inout) :: self !! Encounter list
class(swiftest_encounter), intent(in) :: source !! Source object to copy into
end subroutine symba_util_copy_encounter
end interface

interface util_fill
Expand Down
25 changes: 23 additions & 2 deletions src/symba/symba_util.f90
Original file line number Diff line number Diff line change
Expand Up @@ -188,6 +188,28 @@ module subroutine symba_util_append_tp(self, source, lsource_mask)
end subroutine symba_util_append_tp


module subroutine symba_util_copy_encounter(self, source)
!! author: David A. Minton
!!
!! Copies elements from the source encounter list into self.
implicit none
! Arguments
class(symba_encounter), intent(inout) :: self !! Encounter list
class(swiftest_encounter), intent(in) :: source !! Source object to copy into

select type(source)
class is (symba_encounter)
associate(n => source%nenc)
self%level(1:n) = source%level(1:n)
end associate
end select

call util_copy_encounter(self, source)

return
end subroutine symba_util_copy_encounter


module subroutine symba_util_fill_arr_info(keeps, inserts, lfill_list)
!! author: David A. Minton
!!
Expand Down Expand Up @@ -419,6 +441,7 @@ module subroutine symba_util_rearray_pl(self, system, param)
if (pl_adds%nbody > 0) then
! First store the original plplenc list so we don't remove any of the original encounters
allocate(plplenc_old, source=system%plplenc_list)
call plplenc_old%copy(system%plplenc_list)

! Append the adds to the main pl object
call pl%append(pl_adds, lsource_mask=[(.true., i=1, pl_adds%nbody)])
Expand Down Expand Up @@ -466,8 +489,6 @@ module subroutine symba_util_rearray_pl(self, system, param)
end if
end do
end associate
call move_alloc(plplenc_old, system%plplenc_list)

end if

end associate
Expand Down

0 comments on commit b6d35fb

Please sign in to comment.