Returns the state of the from
frame center w.r.t. the to
frame center,
at the specified ephemeris time et
.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(reference_frame), | intent(in) | :: | from | |||
class(reference_frame), | intent(in) | :: | to | |||
class(ephemeris_class), | intent(inout) | :: | eph |
for ephemeris computations (assumed to have already been initialized) |
||
real(kind=wp), | intent(in) | :: | et |
ephemeris time [sec] |
||
real(kind=wp), | intent(out), | dimension(3) | :: | rc21 |
position of |
|
real(kind=wp), | intent(out), | dimension(3) | :: | vc21 |
velocity of |
|
logical, | intent(out) | :: | status_ok |
true if there were no errors |
subroutine rvcto_rvcfrom_icrf(from, to, eph, et, rc21, vc21, status_ok) !! Returns the state of the `from` frame center w.r.t. the `to` frame center, !! at the specified ephemeris time `et`. implicit none class(reference_frame),intent(in) :: from class(reference_frame),intent(in) :: to class(ephemeris_class),intent(inout) :: eph !! for ephemeris computations (assumed to have already been initialized) real(wp),intent(in) :: et !! ephemeris time [sec] real(wp),dimension(3),intent(out) :: rc21 !! position of `from` frame center w.r.t. `to` frame center real(wp),dimension(3),intent(out) :: vc21 !! velocity of `from` frame center w.r.t. `to` frame center logical,intent(out) :: status_ok !! true if there were no errors ! local variables real(wp),dimension(6) :: rvc1 !! inertial state of `from` frame center w.r.t. `from` primary body real(wp),dimension(6) :: rvc2 !! inertial state of `to` frame center w.r.t. `to` primary body real(wp),dimension(6) :: rvb21 !! inertial state of `from` primary body w.r.t. `to` primary body real(wp),dimension(6) :: rvc21 !! inertial state of `from` frame center w.r.t. `to` frame center ! get TO primary body -> FROM primary body state [inertial] call eph%get_rv(et,from%primary_body,to%primary_body,rvb21,status_ok) if (status_ok) then ! currently, only the two-body rotating frames may be ! centered somewhere other than the primary body. select type (from) class is (two_body_rotating_frame) call from%from_primary_to_center(eph,et,rvc1,status_ok) class default rvc1 = zero end select if (status_ok) then select type (to) class is (two_body_rotating_frame) call to%from_primary_to_center(eph,et,rvc2,status_ok) class default rvc2 = zero end select if (status_ok) then rvc21 = rvb21 + rvc1 - rvc2 rc21 = rvc21(1:3) vc21 = rvc21(4:6) return else write(error_unit,'(A)') 'Error in rvcto_rvcfrom_icrf: '//& 'Could not compute center of TO frame.' end if else write(error_unit,'(A)') 'Error in rvcto_rvcfrom_icrf: '//& 'Could not compute center of FROM frame.' end if else !error write(error_unit,'(A)') 'Error in rvcto_rvcfrom_icrf: '//& 'Could not compute translation.' end if ! we end up here if there was an error: rc21 = zero vc21 = zero end subroutine rvcto_rvcfrom_icrf