from_ijk_to_frame_rv Subroutine

private subroutine from_ijk_to_frame_rv(mu, from_ijk_to_frame, rt_ijk, vt_ijk, r_ijk, v_ijk, dr_frame, dv_frame)

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

Arguments

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

gravitational parameter [km^3/s^2]

procedure(frame_transform_func) :: from_ijk_to_frame

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) :: r_ijk

Chaser IJK absolute position vector [km]

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

Chaser IJK absolute position vector [km]

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

Chaser frame position vector relative to target [km]

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

Chaser frame position vector relative to target [km]


Called by

proc~~from_ijk_to_frame_rv~~CalledByGraph proc~from_ijk_to_frame_rv from_ijk_to_frame_rv proc~from_ijk_to_lvlh_rv from_ijk_to_lvlh_rv proc~from_ijk_to_lvlh_rv->proc~from_ijk_to_frame_rv proc~from_ijk_to_rsw_rv from_ijk_to_rsw_rv proc~from_ijk_to_rsw_rv->proc~from_ijk_to_frame_rv proc~from_ijk_to_vuw_rv from_ijk_to_vuw_rv proc~from_ijk_to_vuw_rv->proc~from_ijk_to_frame_rv interface~from_ijk_to_lvlh from_ijk_to_lvlh interface~from_ijk_to_lvlh->proc~from_ijk_to_lvlh_rv interface~from_ijk_to_rsw from_ijk_to_rsw interface~from_ijk_to_rsw->proc~from_ijk_to_rsw_rv interface~from_ijk_to_vuw from_ijk_to_vuw interface~from_ijk_to_vuw->proc~from_ijk_to_vuw_rv proc~from_lvlh_to_ijk_mat from_lvlh_to_ijk_mat proc~from_lvlh_to_ijk_mat->interface~from_ijk_to_lvlh proc~from_rsw_to_ijk_mat from_rsw_to_ijk_mat proc~from_rsw_to_ijk_mat->interface~from_ijk_to_rsw proc~from_vuw_to_ijk_mat from_vuw_to_ijk_mat proc~from_vuw_to_ijk_mat->interface~from_ijk_to_vuw proc~relative_motion_test relative_motion_test proc~relative_motion_test->interface~from_ijk_to_lvlh proc~relative_motion_test->interface~from_ijk_to_rsw interface~from_lvlh_to_ijk from_lvlh_to_ijk interface~from_lvlh_to_ijk->proc~from_lvlh_to_ijk_mat interface~from_rsw_to_ijk from_rsw_to_ijk interface~from_rsw_to_ijk->proc~from_rsw_to_ijk_mat interface~from_vuw_to_ijk from_vuw_to_ijk interface~from_vuw_to_ijk->proc~from_vuw_to_ijk_mat

Source Code

    subroutine from_ijk_to_frame_rv(mu,from_ijk_to_frame,rt_ijk,vt_ijk,r_ijk,v_ijk,dr_frame,dv_frame)

    implicit none

    real(wp),intent(in)                        :: mu       !! gravitational parameter [km^3/s^2]
    procedure(frame_transform_func)            :: from_ijk_to_frame !! 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)           :: r_ijk    !! Chaser IJK absolute position vector [km]
    real(wp),dimension(3),intent(in)           :: v_ijk    !! Chaser IJK absolute position vector [km]
    real(wp),dimension(3),intent(out)          :: dr_frame !! Chaser frame position vector relative to target [km]
    real(wp),dimension(3),intent(out),optional :: dv_frame !! Chaser frame position vector relative to target [km]

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

    !IJK state of chaser relative to target:
    dr_ijk = r_ijk - rt_ijk    ! [target + delta = chaser]

    if (present(dv_frame)) then

        dv_ijk = v_ijk - vt_ijk ! [target + delta = chaser]

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

        dr_frame = matmul( c, dr_ijk )
        dv_frame = matmul( cdot, dr_ijk ) + matmul( c, dv_ijk )

    else

        call from_ijk_to_frame(mu,r_ijk,v_ijk,c=c)

        dr_frame = matmul( c, dr_ijk )

    end if

    end subroutine from_ijk_to_frame_rv