normalize_variables Subroutine

public subroutine normalize_variables(mu1, mu2, d12, x, m, t, x_crtbp, m_crtbp, t_crtbp)

Convert state in km, km/s units to normalized CRTBP state.

Arguments

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

grav. param. of body 1

real(kind=wp), intent(in) :: mu2

grav. param. of body 2

real(kind=wp), intent(in) :: d12

distance between body 1 and body 2 [km]

real(kind=wp), intent(in), optional, dimension(6) :: x

unnormalized state w.r.t. barycenter [km,km/s]

real(kind=wp), intent(in), optional :: m

unnormalized mass [kg]

real(kind=wp), intent(in), optional :: t

unnormalized time [sec]

real(kind=wp), intent(out), optional, dimension(6) :: x_crtbp

CRTBP normalized state

real(kind=wp), intent(out), optional :: m_crtbp

CRTBP normalized mass

real(kind=wp), intent(out), optional :: t_crtbp

CRTBP normalized time


Source Code

    subroutine normalize_variables(mu1,mu2,d12,x,m,t,x_crtbp,m_crtbp,t_crtbp)

    implicit none

    real(wp),intent(in)                        :: mu1      !! grav. param. of body 1 \( [\mathrm{km}^3/\mathrm{s}^2] \)
    real(wp),intent(in)                        :: mu2      !! grav. param. of body 2 \( [\mathrm{km}^3/\mathrm{s}^2] \)
    real(wp),intent(in)                        :: d12      !! distance between body 1 and body 2 [km]
    real(wp),dimension(6),intent(in),optional  :: x        !! unnormalized state w.r.t. barycenter [km,km/s]
    real(wp),intent(in),optional               :: m        !! unnormalized mass [kg]
    real(wp),intent(in),optional               :: t        !! unnormalized time [sec]
    real(wp),dimension(6),intent(out),optional :: x_crtbp  !! CRTBP normalized state
    real(wp),intent(out),optional              :: m_crtbp  !! CRTBP normalized mass
    real(wp),intent(out),optional              :: t_crtbp  !! CRTBP normalized time

    real(wp) :: tp,n,tu,du,mu

    n = sqrt((mu1+mu2)/d12**3)  !mean motion (rad/sec)
    tu = one/n  !time unit
    du = d12    !distance unit

    if (present(x) .and. present(x_crtbp)) then
        x_crtbp(1:3) = x(1:3) / du       !scale distance
        x_crtbp(4:6) = x(4:6) / (du/tu)  !scale velocity
    end if

    if (present(m) .and. present(m_crtbp)) then
        mu = (mu1 + mu2)/universal_grav_constant    !mass unit
        m_crtbp = m / mu  !scale mass
    end if

    if (present(t) .and. present(t_crtbp)) t_crtbp = t / tu  !scale time

    end subroutine normalize_variables