Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(jpl_ephemeris_splined), | intent(inout) | :: | me | |||
real(kind=wp), | intent(in) | :: | et |
ephemeris time [sec] |
||
type(celestial_body), | intent(in) | :: | targ |
target body |
||
type(celestial_body), | intent(in) | :: | obs |
observer body |
||
real(kind=wp), | intent(out), | dimension(3) | :: | r |
r of targ w.r.t. obs [km] in ICRF frame |
|
logical, | intent(out) | :: | status_ok |
true if there were no problems |
subroutine get_r_splined(me,et,targ,obs,r,status_ok) class(jpl_ephemeris_splined),intent(inout) :: me real(wp),intent(in) :: et !! ephemeris time [sec] type(celestial_body),intent(in) :: targ !! target body type(celestial_body),intent(in) :: obs !! observer body real(wp),dimension(3),intent(out) :: r !! r of targ w.r.t. obs [km] in ICRF frame logical,intent(out) :: status_ok !! true if there were no problems status_ok = .true. if (targ==body_earth .and. obs==body_moon) then r = me%earth_eph_interface%get_r(et) elseif (targ==body_sun .and. obs==body_moon) then r = me%sun_eph_interface%get_r(et) elseif (targ==body_ssb .and. obs==body_moon) then r = me%ssb_eph_interface%get_r(et) elseif (targ==body_moon .and. obs==body_earth) then ! inverse are negative r = -me%earth_eph_interface%get_r(et) elseif (targ==body_moon .and. obs==body_sun) then r = -me%sun_eph_interface%get_r(et) elseif (targ==body_moon .and. obs==body_ssb) then r = -me%ssb_eph_interface%get_r(et) elseif (targ==body_sun .and. obs==body_ssb) then ! for this one we subtract these ! ssb -> sun = ssb -> moon + moon -> sun r = me%sun_eph_interface%get_r(et) - me%ssb_eph_interface%get_r(et) elseif (targ==body_earth .and. obs==body_ssb) then ! for this one we subtract these ! ssb -> earth = ssb -> moon + moon -> earth r = me%earth_eph_interface%get_r(et) - me%ssb_eph_interface%get_r(et) else write(*,*) 'targ = ', targ write(*,*) 'obs = ', obs error stop 'error in get_r_splined: this combo has not been splined' ! or could call me%jpl_ephemeris%get_rv(et,targ,obs,rv,status_ok); r = rv(1:3) end if end subroutine get_r_splined