uhat_dot Function

public pure function uhat_dot(u, udot) result(uhatd)

Time derivative of a unit vector.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in), dimension(3) :: u

vector [u]

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

derivative of vector [du/dt]

Return Value real(kind=wp), dimension(3)

derivative of unit vector [d(uhat)/dt]


Called by

proc~~uhat_dot~~CalledByGraph proc~uhat_dot vector_module::uhat_dot proc~from_ijk_to_lvlh_mat relative_motion_module::from_ijk_to_lvlh_mat proc~from_ijk_to_lvlh_mat->proc~uhat_dot proc~from_ijk_to_rsw_mat relative_motion_module::from_ijk_to_rsw_mat proc~from_ijk_to_rsw_mat->proc~uhat_dot interface~from_ijk_to_lvlh relative_motion_module::from_ijk_to_lvlh 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 interface~from_ijk_to_rsw relative_motion_module::from_ijk_to_rsw interface~from_ijk_to_rsw->proc~from_ijk_to_rsw_mat proc~from_ijk_to_rsw_rv relative_motion_module::from_ijk_to_rsw_rv interface~from_ijk_to_rsw->proc~from_ijk_to_rsw_rv proc~from_ijk_to_lvlh_rv->interface~from_ijk_to_lvlh proc~from_ijk_to_rsw_rv->interface~from_ijk_to_rsw proc~from_lvlh_to_ijk_mat relative_motion_module::from_lvlh_to_ijk_mat proc~from_lvlh_to_ijk_mat->interface~from_ijk_to_lvlh proc~from_rsw_to_ijk_mat relative_motion_module::from_rsw_to_ijk_mat proc~from_rsw_to_ijk_mat->interface~from_ijk_to_rsw proc~relative_motion_test relative_motion_module::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 relative_motion_module::from_lvlh_to_ijk interface~from_lvlh_to_ijk->proc~from_lvlh_to_ijk_mat proc~from_lvlh_to_ijk_rv relative_motion_module::from_lvlh_to_ijk_rv interface~from_lvlh_to_ijk->proc~from_lvlh_to_ijk_rv interface~from_rsw_to_ijk relative_motion_module::from_rsw_to_ijk interface~from_rsw_to_ijk->proc~from_rsw_to_ijk_mat proc~from_rsw_to_ijk_rv relative_motion_module::from_rsw_to_ijk_rv interface~from_rsw_to_ijk->proc~from_rsw_to_ijk_rv proc~from_lvlh_to_ijk_rv->interface~from_lvlh_to_ijk proc~from_rsw_to_ijk_rv->interface~from_rsw_to_ijk

Source Code

    pure function uhat_dot(u,udot) result(uhatd)

    implicit none

    real(wp),dimension(3),intent(in) :: u      !! vector [`u`]
    real(wp),dimension(3),intent(in) :: udot   !! derivative of vector [`du/dt`]
    real(wp),dimension(3)            :: uhatd  !! derivative of unit vector [`d(uhat)/dt`]

    real(wp)              :: umag  !! vector magnitude
    real(wp),dimension(3) :: uhat  !! unit vector

    umag = norm2(u)

    if (umag == zero) then  !singularity
        uhatd = zero
    else
        uhat = u / umag
        uhatd = ( udot - dot_product(uhat,udot)*uhat ) / umag
    end if

    end function uhat_dot