periapsis_raise_maneuver Function

private pure function periapsis_raise_maneuver(rv, target_rp) result(dv)

Compute the maneuver at apoapsis to raise periapsis to the specified value.

Arguments

TypeIntentOptionalAttributesName
real(kind=wp), intent(in), dimension(6):: rv

[position,velocity] vector

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

the rp value to target

Return Value real(kind=wp)

the maneuver to perform at apoapsis to target target_rp


Calls

proc~~periapsis_raise_maneuver~~CallsGraph proc~periapsis_raise_maneuver periapsis_raise_maneuver rv_to_orbital_elements rv_to_orbital_elements proc~periapsis_raise_maneuver->rv_to_orbital_elements periapsis_apoapsis periapsis_apoapsis proc~periapsis_raise_maneuver->periapsis_apoapsis

Called by

proc~~periapsis_raise_maneuver~~CalledByGraph proc~periapsis_raise_maneuver periapsis_raise_maneuver proc~altitude_maintenance altitude_maintenance proc~altitude_maintenance->proc~periapsis_raise_maneuver

Contents


Source Code

    pure function periapsis_raise_maneuver(rv,target_rp) result(dv)

    implicit none

    real(wp),dimension(6),intent(in) :: rv !! [position,velocity] vector
    real(wp),intent(in) :: target_rp  !! the rp value to target
    real(wp) :: dv !! the maneuver to perform at apoapsis to target `target_rp`

    real(wp),dimension(3) :: r !! position vector
    real(wp),dimension(3) :: v !! velocity vector
    real(wp) :: a,p,ecc,inc,raan,aop,tru !! orbital elements
    real(wp) :: rp1,ra1,vp1,va1,va2  !! periapsis/apoapsis pos/vel magnitudes

    r = rv(1:3)
    v = rv(4:6)

    call rv_to_orbital_elements(body_moon%mu,r,v,p,ecc,inc,raan,aop,tru)
    a = p / (one - ecc*ecc) ! compute semi-major axis
    call periapsis_apoapsis(body_moon%mu,a,ecc,rp1,ra1,vp1,va1)

    ! apoapsis velocity for periapsis radius of target_rp
    va2 = sqrt( two * body_moon%mu * target_rp/(ra1*(target_rp+ra1)) )

    ! delta-v to raise periapsis back to target rp
    dv = va2 - va1

    end function periapsis_raise_maneuver