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

Commit

Permalink
Added in xh2xb coordinate transform method
Browse files Browse the repository at this point in the history
  • Loading branch information
daminton committed Aug 16, 2021
1 parent 412c1dd commit 3b7603f
Show file tree
Hide file tree
Showing 3 changed files with 44 additions and 1 deletion.
7 changes: 7 additions & 0 deletions src/modules/swiftest_classes.f90
Original file line number Diff line number Diff line change
Expand Up @@ -226,6 +226,7 @@ module swiftest_classes
procedure :: append => util_append_pl !! Appends elements from one structure to another
procedure :: h2b => util_coord_h2b_pl !! Convert massive bodies from heliocentric to barycentric coordinates (position and velocity)
procedure :: b2h => util_coord_b2h_pl !! Convert massive bodies from barycentric to heliocentric coordinates (position and velocity)
procedure :: xh2xb => util_coord_xh2xb_pl !! Convert position vectors of massive bodies from heliocentric to barycentric coordinates (position and velocity)
procedure :: fill => util_fill_pl !! "Fills" bodies from one object into another depending on the results of a mask (uses the UNPACK intrinsic)
procedure :: resize => util_resize_pl !! Checks the current size of a Swiftest body against the requested size and resizes it if it is too small.
procedure :: set_beg_end => util_set_beg_end_pl !! Sets the beginning and ending positions and velocities of planets.
Expand Down Expand Up @@ -936,6 +937,12 @@ module subroutine util_coord_h2b_tp(self, cb)
class(swiftest_cb), intent(in) :: cb !! Swiftest central body object
end subroutine util_coord_h2b_tp

module subroutine util_coord_xh2xb_pl(self, cb)
implicit none
class(swiftest_pl), intent(inout) :: self !! Swiftest massive body object
class(swiftest_cb), intent(inout) :: cb !! Swiftest central body object
end subroutine util_coord_xh2xb_pl

module subroutine util_copy_encounter(self, source)
implicit none
class(swiftest_encounter), intent(inout) :: self !! Encounter list
Expand Down
2 changes: 1 addition & 1 deletion src/symba/symba_collision.f90
Original file line number Diff line number Diff line change
Expand Up @@ -1076,7 +1076,7 @@ module subroutine symba_collision_resolve_plplenc(self, system, param, t)

write(*, *) "Collision between massive bodies detected at time t = ", t
call pl%vb2vh(system%cb)
call pl%h2b(system%cb)
call pl%xh2xb(system%cb)
if (param%lfragmentation) then
call plplcollision_list%resolve_fragmentations(system, param)
else
Expand Down
36 changes: 36 additions & 0 deletions src/util/util_coord.f90
Original file line number Diff line number Diff line change
Expand Up @@ -42,6 +42,42 @@ module subroutine util_coord_h2b_pl(self, cb)
end subroutine util_coord_h2b_pl


module subroutine util_coord_xh2xb_pl(self, cb)
!! author: David A. Minton
!!
!! Convert position vectors of massive bodies from heliocentric to barycentric coordinates (position and velocity)
!!
!! Adapted from David E. Kaufmann's Swifter routine coord_h2b.f90
!! Adapted from Hal Levison's Swift routine coord_h2b.f
implicit none
! Arguments
class(swiftest_pl), intent(inout) :: self !! Swiftest massive body object
class(swiftest_cb), intent(inout) :: cb !! Swiftest central body object
! Internals
integer(I4B) :: i
real(DP) :: Gmtot
real(DP), dimension(NDIM) :: xtmp

if (self%nbody == 0) return
associate(pl => self, npl => self%nbody)
Gmtot = cb%Gmass
xtmp(:) = 0.0_DP
do i = 1, npl
if (pl%status(i) == INACTIVE) cycle
Gmtot = Gmtot + pl%Gmass(i)
xtmp(:) = xtmp(:) + pl%Gmass(i) * pl%xh(:,i)
end do
cb%xb(:) = -xtmp(:) / Gmtot
do i = 1, npl
if (pl%status(i) == INACTIVE) cycle
pl%xb(:,i) = pl%xh(:,i) + cb%xb(:)
end do
end associate

return
end subroutine util_coord_xh2xb_pl


module subroutine util_coord_h2b_tp(self, cb)
!! author: David A. Minton
!!
Expand Down

0 comments on commit 3b7603f

Please sign in to comment.