Convert modified equinoctial elements (posigrade formulation) to Cartesian coordinates.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=wp), | intent(in) | :: | mu |
central body gravitational parameter () |
||
real(kind=wp), | intent(in), | dimension(6) | :: | evec |
Modified equinoctial element vector |
|
real(kind=wp), | intent(out), | dimension(6) | :: | rv |
Cartesian state vector |
subroutine equinoctial_to_cartesian(mu,evec,rv) implicit none real(wp),intent(in) :: mu !! central body gravitational parameter (\(km^3/s^2\)) real(wp),dimension(6),intent(in) :: evec !! Modified equinoctial element vector real(wp),dimension(6),intent(out) :: rv !! Cartesian state vector real(wp) :: p,f,g,h,k,L,s2,r,w,cL,sL,smp,hh,kk,tkh real(wp) :: x,y,xdot,ydot real(wp),dimension(3) :: fhat,ghat p = evec(1) f = evec(2) g = evec(3) h = evec(4) k = evec(5) L = evec(6) kk = k*k hh = h*h tkh = two*k*h s2 = one + hh + kk cL = cos(L) sL = sin(L) w = one + f*cL + g*sL r = p/w smp = sqrt(mu/p) fhat(1) = one-kk+hh fhat(2) = tkh fhat(3) = -two*k ghat(1) = tkh ghat(2) = one+kk-hh ghat(3) = two*h fhat = fhat/s2 ghat = ghat/s2 x = r*cL y = r*sL xdot = -smp * (g + sL) ydot = smp * (f + cL) rv(1:3) = x*fhat + y*ghat rv(4:6) = xdot*fhat + ydot*ghat end subroutine equinoctial_to_cartesian