Return the state of the targ
body relative to
the obs
body, in the inertial frame [ICRF].
This is the function that can be used in the ephemeris_class.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(standish_ephemeris), | 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(6) | :: | rv |
state of targ w.r.t. obs [km, km/s] |
|
logical, | intent(out) | :: | status_ok |
true if there were no problems |
subroutine standish_rv_func(me,et,targ,obs,rv,status_ok) implicit none class(standish_ephemeris),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(6),intent(out) :: rv !! state of targ w.r.t. obs [km, km/s] logical,intent(out) :: status_ok !! true if there were no problems real(wp) :: jd integer :: itbl_targ, itbl_obs integer :: itarg integer :: iobs real(wp), dimension(6) :: targ_rv_au real(wp), dimension(6) :: obs_rv_au real(wp),dimension(6) :: rv_ecliptic integer,parameter :: naif_id_sun = 10 !! NAIF ID for the Sun if (targ==obs) then rv = zero return end if jd = et_to_jd(et) if (targ%id==naif_id_sun) then itarg = -1 ! dummy else itarg = spice_id_to_standish_id(targ%id) end if if (obs%id==naif_id_sun) then iobs = -1 ! dummy else iobs = spice_id_to_standish_id(obs%id) end if if (itarg==0 .or. iobs==0) then ! if the bodies were not found in this ephemeris rv = zero status_ok = .false. else if (targ%id/=naif_id_sun) then ! targ w.r.t sun [j2000 ecliptic] call helio (itarg, jd, targ_rv_au, itbl_targ) else targ_rv_au = zero itbl_targ = 3 ! dummy value end if if (obs%id/=naif_id_sun) then ! obs w.r.t sun [j2000 ecliptic] call helio (iobs, jd, obs_rv_au, itbl_obs ) else obs_rv_au = zero itbl_obs = 3 ! dummy value end if if (itbl_targ>0 .and. itbl_obs>0) then ! vector from obs to targ [j2000 ecliptic] rv_ecliptic = targ_rv_au - obs_rv_au ! convert to ICRF: call ec2eq (rv_ecliptic, rv) ! Convert to km, km/s: rv = rv * au2m * m2km rv(4:6) = rv(4:6) / (year2day*day2sec) status_ok = .true. else ! out of range of the ephemeris: rv = zero status_ok = .false. end if end if end subroutine standish_rv_func