from_frame_to_ijk_rv Subroutine

private subroutine from_frame_to_ijk_rv(mu, from_frame_to_ijk, rt_ijk, vt_ijk, dr_frame, dv_frame, r_ijk, v_ijk)

Transform a position (and optionally velocity) vector from a specified relative frame to IJK.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: mu

gravitational parameter [km^3/s^2]

procedure(frame_transform_func) :: from_frame_to_ijk

function to compute the transformation matrices

real(kind=wp), intent(in), dimension(3) :: rt_ijk

Target IJK absolute position vector [km]

real(kind=wp), intent(in), dimension(3) :: vt_ijk

Target IJK absolute position vector [km]

real(kind=wp), intent(in), dimension(3) :: dr_frame

Chaser frame position vector relative to target [km]

real(kind=wp), intent(in), dimension(3) :: dv_frame

Chaser frame position vector relative to target [km]

real(kind=wp), intent(out), dimension(3) :: r_ijk

Chaser IJK absolute position vector [km]

real(kind=wp), intent(out), optional, dimension(3) :: v_ijk

Chaser IJK absolute position vector [km]


Called by

proc~~from_frame_to_ijk_rv~~CalledByGraph proc~from_frame_to_ijk_rv from_frame_to_ijk_rv proc~from_lvlh_to_ijk_rv from_lvlh_to_ijk_rv proc~from_lvlh_to_ijk_rv->proc~from_frame_to_ijk_rv proc~from_rsw_to_ijk_rv from_rsw_to_ijk_rv proc~from_rsw_to_ijk_rv->proc~from_frame_to_ijk_rv proc~from_vuw_to_ijk_rv from_vuw_to_ijk_rv proc~from_vuw_to_ijk_rv->proc~from_frame_to_ijk_rv interface~from_lvlh_to_ijk from_lvlh_to_ijk interface~from_lvlh_to_ijk->proc~from_lvlh_to_ijk_rv interface~from_rsw_to_ijk from_rsw_to_ijk interface~from_rsw_to_ijk->proc~from_rsw_to_ijk_rv interface~from_vuw_to_ijk from_vuw_to_ijk interface~from_vuw_to_ijk->proc~from_vuw_to_ijk_rv

Source Code

    subroutine from_frame_to_ijk_rv(mu,from_frame_to_ijk,rt_ijk,vt_ijk,dr_frame,dv_frame,r_ijk,v_ijk)

    implicit none

    real(wp),intent(in)                         :: mu       !! gravitational parameter [km^3/s^2]
    procedure(frame_transform_func)             :: from_frame_to_ijk !! function to compute the transformation matrices
    real(wp),dimension(3),intent(in)            :: rt_ijk    !! Target IJK absolute position vector [km]
    real(wp),dimension(3),intent(in)            :: vt_ijk    !! Target IJK absolute position vector [km]
    real(wp),dimension(3),intent(in)            :: dr_frame  !! Chaser frame position vector relative to target [km]
    real(wp),dimension(3),intent(in)            :: dv_frame  !! Chaser frame position vector relative to target [km]
    real(wp),dimension(3),intent(out)           :: r_ijk     !! Chaser IJK absolute position vector [km]
    real(wp),dimension(3),intent(out),optional  :: v_ijk     !! Chaser IJK absolute position vector [km]

    real(wp),dimension(3,3) :: c
    real(wp),dimension(3,3) :: cdot

    if (present(v_ijk)) then

        call from_frame_to_ijk(mu,rt_ijk,vt_ijk,c=c,cdot=cdot)

        !chaser = target + delta:
        r_ijk = rt_ijk + matmul( c, dr_frame )
        v_ijk = vt_ijk + matmul( cdot, dr_frame ) + matmul( c, dv_frame )

    else

        call from_frame_to_ijk(mu,rt_ijk,vt_ijk,c=c)

        r_ijk = rt_ijk + matmul( c, dr_frame )

    end if

    end subroutine from_frame_to_ijk_rv