from_lvlh_to_ijk_rv Subroutine

private subroutine from_lvlh_to_ijk_rv(rt_ijk, vt_ijk, dr_lvlh, dv_lvlh, r_ijk, v_ijk)

Transform a position (and optionally velocity) vector from LVLH to IJK.

Arguments

Type IntentOptional Attributes Name
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_lvlh

Chaser LVLH position vector relative to target [km]

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

Chaser LVLH 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]


Calls

proc~~from_lvlh_to_ijk_rv~~CallsGraph proc~from_lvlh_to_ijk_rv relative_motion_module::from_lvlh_to_ijk_rv interface~from_lvlh_to_ijk relative_motion_module::from_lvlh_to_ijk proc~from_lvlh_to_ijk_rv->interface~from_lvlh_to_ijk interface~from_lvlh_to_ijk->proc~from_lvlh_to_ijk_rv proc~from_lvlh_to_ijk_mat relative_motion_module::from_lvlh_to_ijk_mat interface~from_lvlh_to_ijk->proc~from_lvlh_to_ijk_mat interface~from_ijk_to_lvlh relative_motion_module::from_ijk_to_lvlh proc~from_lvlh_to_ijk_mat->interface~from_ijk_to_lvlh proc~from_ijk_to_lvlh_mat relative_motion_module::from_ijk_to_lvlh_mat interface~from_ijk_to_lvlh->proc~from_ijk_to_lvlh_mat proc~from_ijk_to_lvlh_rv relative_motion_module::from_ijk_to_lvlh_rv interface~from_ijk_to_lvlh->proc~from_ijk_to_lvlh_rv proc~cross vector_module::cross proc~from_ijk_to_lvlh_mat->proc~cross proc~uhat_dot vector_module::uhat_dot proc~from_ijk_to_lvlh_mat->proc~uhat_dot proc~unit vector_module::unit proc~from_ijk_to_lvlh_mat->proc~unit proc~from_ijk_to_lvlh_rv->interface~from_ijk_to_lvlh

Called by

proc~~from_lvlh_to_ijk_rv~~CalledByGraph proc~from_lvlh_to_ijk_rv relative_motion_module::from_lvlh_to_ijk_rv interface~from_lvlh_to_ijk relative_motion_module::from_lvlh_to_ijk proc~from_lvlh_to_ijk_rv->interface~from_lvlh_to_ijk interface~from_lvlh_to_ijk->proc~from_lvlh_to_ijk_rv

Source Code

    subroutine from_lvlh_to_ijk_rv(rt_ijk,vt_ijk,dr_lvlh,dv_lvlh,r_ijk,v_ijk)

    implicit none

    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_lvlh  !! Chaser LVLH position vector relative to target [km]
    real(wp),dimension(3),intent(in)            :: dv_lvlh  !! Chaser LVLH 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_lvlh_to_ijk(rt_ijk,vt_ijk,c=c,cdot=cdot)

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

    else

        call from_lvlh_to_ijk(rt_ijk,vt_ijk,c=c)

        r_ijk = rt_ijk + matmul( c, dr_lvlh )

    end if

    end subroutine from_lvlh_to_ijk_rv