def gooding_coe(k, p, ecc, inc, raan, argp, nu, tof, numiter=150, rtol=1e-8): # TODO: parabolic and hyperbolic not implemented cases if ecc >= 1.0: raise NotImplementedError( "Parabolic/Hyperbolic cases still not implemented in gooding.") M0 = E_to_M(nu_to_E(nu, ecc), ecc) semi_axis_a = p / (1 - ecc**2) n = np.sqrt(k / np.abs(semi_axis_a)**3) M = M0 + n * tof # Start the computation n = 0 c = ecc * np.cos(M) s = ecc * np.sin(M) psi = s / np.sqrt(1 - 2 * c + ecc**2) f = 1.0 while f**2 >= rtol and n <= numiter: xi = np.cos(psi) eta = np.sin(psi) fd = (1 - c * xi) + s * eta fdd = c * eta + s * xi f = psi - fdd psi = psi - f * fd / (fd**2 - 0.5 * f * fdd) n += 1 E = M + psi return E_to_nu(E, ecc)
def markley_coe(k, p, ecc, inc, raan, argp, nu, tof): M0 = E_to_M(nu_to_E(nu, ecc), ecc) a = p / (1 - ecc**2) n = np.sqrt(k / a**3) M = M0 + n * tof # Range between -pi and pi M = (M + np.pi) % (2 * np.pi) - np.pi # Equation (20) alpha = (3 * np.pi**2 + 1.6 * (np.pi - np.abs(M)) / (1 + ecc)) / (np.pi**2 - 6) # Equation (5) d = 3 * (1 - ecc) + alpha * ecc # Equation (9) q = 2 * alpha * d * (1 - ecc) - M**2 # Equation (10) r = 3 * alpha * d * (d - 1 + ecc) * M + M**3 # Equation (14) w = (np.abs(r) + np.sqrt(q**3 + r**2))**(2 / 3) # Equation (15) E = (2 * r * w / (w**2 + w * q + q**2) + M) / d # Equation (26) f0 = _kepler_equation(E, M, ecc) f1 = _kepler_equation_prime(E, M, ecc) f2 = ecc * np.sin(E) f3 = ecc * np.cos(E) f4 = -f2 # Equation (22) delta3 = -f0 / (f1 - 0.5 * f0 * f2 / f1) delta4 = -f0 / (f1 + 0.5 * delta3 * f2 + 1 / 6 * delta3**2 * f3) delta5 = -f0 / (f1 + 0.5 * delta4 * f2 + 1 / 6 * delta4**2 * f3 + 1 / 24 * delta4**3 * f4) E += delta5 nu = E_to_nu(E, ecc) return nu
def pimienta_coe(k, p, ecc, inc, raan, argp, nu, tof): q = p / (1 + ecc) # TODO: Do something to increase parabolic accuracy? n = np.sqrt(k * (1 - ecc)**3 / q**3) M0 = E_to_M(nu_to_E(nu, ecc), ecc) M = M0 + n * tof # Equation (32a), (32b), (32c) and (32d) c3 = 5 / 2 + 560 * ecc a = 15 * (1 - ecc) / c3 b = -M / c3 y = np.sqrt(b**2 / 4 + a**3 / 27) # Equation (33) x_bar = (-b / 2 + y)**(1 / 3) - (b / 2 + y)**(1 / 3) # Coefficients from equations (34a) and (34b) c15 = 3003 / 14336 + 16384 * ecc c13 = 3465 / 13312 - 61440 * ecc c11 = 945 / 2816 + 92160 * ecc c9 = 175 / 384 - 70400 * ecc c7 = 75 / 112 + 28800 * ecc c5 = 9 / 8 - 6048 * ecc # Precompute x_bar powers, equations (35a) to (35d) x_bar2 = x_bar**2 x_bar3 = x_bar2 * x_bar x_bar4 = x_bar3 * x_bar x_bar5 = x_bar4 * x_bar x_bar6 = x_bar5 * x_bar x_bar7 = x_bar6 * x_bar x_bar8 = x_bar7 * x_bar x_bar9 = x_bar8 * x_bar x_bar10 = x_bar9 * x_bar x_bar11 = x_bar10 * x_bar x_bar12 = x_bar11 * x_bar x_bar13 = x_bar12 * x_bar x_bar14 = x_bar13 * x_bar x_bar15 = x_bar14 * x_bar # Function f and its derivatives are given by all the (36) equation set f = (c15 * x_bar15 + c13 * x_bar13 + c11 * x_bar11 + c9 * x_bar9 + c7 * x_bar7 + c5 * x_bar5 + c3 * x_bar3 + 15 * (1 - ecc) * x_bar - M) f1 = (15 * c15 * x_bar14 + 13 * c13 * x_bar12 + 11 * c11 * x_bar10 + 9 * c9 * x_bar8 + 7 * c7 * x_bar6 + 5 * c5 * x_bar4 + 3 * c3 * x_bar2 + 15 * (1 - ecc)) f2 = (210 * c15 * x_bar13 + 156 * c13 * x_bar11 + 110 * c11 * x_bar9 + 72 * c9 * x_bar7 + 42 * c7 * x_bar5 + 20 * c5 * x_bar3 + 6 * c3 * x_bar) f3 = (2730 * c15 * x_bar12 + 1716 * c13 * x_bar10 + 990 * c11 * x_bar8 + 504 * c9 * x_bar6 + 210 * c7 * x_bar4 + 60 * c5 * x_bar2 + 6 * c3) f4 = (32760 * c15 * x_bar11 + 17160 * c13 * x_bar9 + 7920 * c11 * x_bar7 + 3024 * c9 * x_bar5 + 840 * c7 * x_bar3 + 120 * c5 * x_bar) f5 = (360360 * c15 * x_bar10 + 154440 * c13 * x_bar8 + 55440 * c11 * x_bar6 + 15120 * c9 * x_bar4 + 2520 * c7 * x_bar2 + 120 * c5) f6 = (3603600 * c15 * x_bar9 + 1235520 * c13 * x_bar7 + 332640 * c11 * x_bar5 + 60480 * c9 * x_bar3 + 5040 * c7 * x_bar) f7 = (32432400 * c15 * x_bar8 + 8648640 * c13 * x_bar6 + 1663200 * c11 * x_bar4 + 181440 * c9 * x_bar2 + 5040 * c7) f8 = (259459200 * c15 * x_bar7 + 51891840 * c13 * x_bar5 + 6652800 * c11 * x_bar3 + 362880 * c9 * x_bar) f9 = (1.8162144e9 * c15 * x_bar6 + 259459200 * c13 * x_bar4 + 19958400 * c11 * x_bar2 + 362880 * c9) f10 = (1.08972864e10 * c15 * x_bar5 + 1.0378368e9 * c13 * x_bar3 + 39916800 * c11 * x_bar) f11 = 5.4486432e10 * c15 * x_bar4 + 3.1135104e9 * c13 * x_bar2 + 39916800 * c11 f12 = 2.17945728e11 * c15 * x_bar3 + 6.2270208e9 * c13 * x_bar f13 = 6.53837184 * c15 * x_bar2 + 6.2270208e9 * c13 f14 = 1.307674368e13 * c15 * x_bar f15 = 1.307674368e13 * c15 # Solving g parameters defined by equations (37a), (37b), (37c) and (37d) g1 = 1 / 2 g2 = 1 / 6 g3 = 1 / 24 g4 = 1 / 120 g5 = 1 / 720 g6 = 1 / 5040 g7 = 1 / 40320 g8 = 1 / 362880 g9 = 1 / 3628800 g10 = 1 / 39916800 g11 = 1 / 479001600 g12 = 1 / 6.2270208e9 g13 = 1 / 8.71782912e10 g14 = 1 / 1.307674368e12 # Solving for the u_{i} and h_{i} variables defined by equation (38) u1 = -f / f1 h2 = f1 + g1 * u1 * f2 u2 = -f / h2 h3 = f1 + g1 * u2 * f2 + g2 * u2**2 * f3 u3 = -f / h3 h4 = f1 + g1 * u3 * f2 + g2 * u3**2 * f3 + g3 * u3**3 * f4 u4 = -f / h4 h5 = f1 + g1 * u4 * f2 + g2 * u4**2 * f3 + g3 * u4**3 * f4 + g4 * u4**4 * f5 u5 = -f / h5 h6 = (f1 + g1 * u5 * f2 + g2 * u5**2 * f3 + g3 * u5**3 * f4 + g4 * u5**4 * f5 + g5 * u5**5 * f6) u6 = -f / h6 h7 = (f1 + g1 * u6 * f2 + g2 * u6**2 * f3 + g3 * u6**3 * f4 + g4 * u6**4 * f5 + g5 * u6**5 * f6 + g6 * u6**6 * f7) u7 = -f / h7 h8 = (f1 + g1 * u7 * f2 + g2 * u7**2 * f3 + g3 * u7**3 * f4 + g4 * u7**4 * f5 + g5 * u7**5 * f6 + g6 * u7**6 * f7 + g7 * u7**7 * f8) u8 = -f / h8 h9 = (f1 + g1 * u8 * f2 + g2 * u8**2 * f3 + g3 * u8**3 * f4 + g4 * u8**4 * f5 + g5 * u8**5 * f6 + g6 * u8**6 * f7 + g7 * u8**7 * f8 + g8 * u8**8 * f9) u9 = -f / h9 h10 = (f1 + g1 * u9 * f2 + g2 * u9**2 * f3 + g3 * u9**3 * f4 + g4 * u9**4 * f5 + g5 * u9**5 * f6 + g6 * u9**6 * f7 + g7 * u9**7 * f8 + g8 * u9**8 * f9 + g9 * u9**9 * f10) u10 = -f / h10 h11 = (f1 + g1 * u10 * f2 + g2 * u10**2 * f3 + g3 * u10**3 * f4 + g4 * u10**4 * f5 + g5 * u10**5 * f6 + g6 * u10**6 * f7 + g7 * u10**7 * f8 + g8 * u10**8 * f9 + g9 * u10**9 * f10 + g10 * u10**10 * f11) u11 = -f / h11 h12 = (f1 + g1 * u11 * f2 + g2 * u11**2 * f3 + g3 * u11**3 * f4 + g4 * u11**4 * f5 + g5 * u11**5 * f6 + g6 * u11**6 * f7 + g7 * u11**7 * f8 + g8 * u11**8 * f9 + g9 * u11**9 * f10 + g10 * u11**10 * f11 + g11 * u11**11 * f12) u12 = -f / h12 h13 = (f1 + g1 * u12 * f2 + g2 * u12**2 * f3 + g3 * u12**3 * f4 + g4 * u12**4 * f5 + g5 * u12**5 * f6 + g6 * u12**6 * f7 + g7 * u12**7 * f8 + g8 * u12**8 * f9 + g9 * u12**9 * f10 + g10 * u12**10 * f11 + g11 * u12**11 * f12 + g12 * u12**12 * f13) u13 = -f / h13 h14 = (f1 + g1 * u13 * f2 + g2 * u13**2 * f3 + g3 * u13**3 * f4 + g4 * u13**4 * f5 + g5 * u13**5 * f6 + g6 * u13**6 * f7 + g7 * u13**7 * f8 + g8 * u13**8 * f9 + g9 * u13**9 * f10 + g10 * u13**10 * f11 + g11 * u13**11 * f12 + g12 * u13**12 * f13 + g13 * u13**13 * f14) u14 = -f / h14 h15 = (f1 + g1 * u14 * f2 + g2 * u14**2 * f3 + g3 * u14**3 * f4 + g4 * u14**4 * f5 + g5 * u14**5 * f6 + g6 * u14**6 * f7 + g7 * u14**7 * f8 + g8 * u14**8 * f9 + g9 * u14**9 * f10 + g10 * u14**10 * f11 + g11 * u14**11 * f12 + g12 * u14**12 * f13 + g13 * u14**13 * f14 + g14 * u14**14 * f15) u15 = -f / h15 # Solving for x x = x_bar + u15 w = x - 0.01171875 * x**17 / (1 + ecc) # Solving for the true anomaly from eccentricity anomaly E = M + ecc * (-16384 * w**15 + 61440 * w**13 - 92160 * w**11 + 70400 * w**9 - 28800 * w**7 + 6048 * w**5 - 560 * w**3 + 15 * w) return E_to_nu(E, ecc)
def recseries_coe(k, p, ecc, inc, raan, argp, nu, tof, method="rtol", order=8, numiter=100, rtol=1e-8): # semi-major axis semi_axis_a = p / (1 - ecc**2) # mean angular motion n = np.sqrt(k / np.abs(semi_axis_a)**3) if ecc == 0: # Solving for circular orbit # compute initial mean anoamly M0 = nu # For circular orbit (M = E = nu) # final mean anaomaly M = M0 + n * tof # snapping anomaly to [0,pi] range nu = M - 2 * np.pi * np.floor(M / 2 / np.pi) return nu elif ecc < 1.0: # Solving for elliptical orbit # compute initial mean anoamly M0 = E_to_M(nu_to_E(nu, ecc), ecc) # final mean anaomaly M = M0 + n * tof # snapping anomaly to [0,pi] range M = M - 2 * np.pi * np.floor(M / 2 / np.pi) # set recursion iteration if method == "rtol": Niter = numiter elif method == "order": Niter = order else: raise ValueError( "Unknown recursion termination method ('rtol','order').") # compute eccentric anomaly through recursive series E = M + ecc # Using initial guess from vallado to improve convergence for i in range(0, Niter): En = M + ecc * np.sin(E) # check for break condition if method == "rtol" and (abs(En - E) / abs(E)) < rtol: break E = En return E_to_nu(E, ecc) else: # Parabolic/Hyperbolic orbits are not supported raise ValueError("Parabolic/Hyperbolic orbits not supported.") return nu
def mikkola_coe(k, p, ecc, inc, raan, argp, nu, tof): a = p / (1 - ecc**2) n = np.sqrt(k / np.abs(a)**3) # Solve for specific geometrical case if ecc < 1.0: # Equation (9a) alpha = (1 - ecc) / (4 * ecc + 1 / 2) M0 = E_to_M(nu_to_E(nu, ecc), ecc) else: alpha = (ecc - 1) / (4 * ecc + 1 / 2) M0 = F_to_M(nu_to_F(nu, ecc), ecc) M = M0 + n * tof beta = M / 2 / (4 * ecc + 1 / 2) # Equation (9b) if beta >= 0: z = (beta + np.sqrt(beta**2 + alpha**3))**(1 / 3) else: z = (beta - np.sqrt(beta**2 + alpha**3))**(1 / 3) s = z - alpha / z # Apply initial correction if ecc < 1.0: ds = -0.078 * s**5 / (1 + ecc) else: ds = 0.071 * s**5 / (1 + 0.45 * s**2) / (1 + 4 * s**2) / ecc s += ds # Solving for the true anomaly if ecc < 1.0: E = M + ecc * (3 * s - 4 * s**3) f = E - ecc * np.sin(E) - M f1 = 1.0 - ecc * np.cos(E) f2 = ecc * np.sin(E) f3 = ecc * np.cos(E) f4 = -f2 f5 = -f3 else: E = 3 * np.log(s + np.sqrt(1 + s**2)) f = -E + ecc * np.sinh(E) - M f1 = -1.0 + ecc * np.cosh(E) f2 = ecc * np.sinh(E) f3 = ecc * np.cosh(E) f4 = f2 f5 = f3 # Apply Taylor expansion u1 = -f / f1 u2 = -f / (f1 + 0.5 * f2 * u1) u3 = -f / (f1 + 0.5 * f2 * u2 + (1.0 / 6.0) * f3 * u2**2) u4 = -f / (f1 + 0.5 * f2 * u3 + (1.0 / 6.0) * f3 * u3**2 + (1.0 / 24.0) * f4 * (u3**3)) u5 = -f / (f1 + f2 * u4 / 2 + f3 * (u4 * u4) / 6.0 + f4 * (u4 * u4 * u4) / 24.0 + f5 * (u4 * u4 * u4 * u4) / 120.0) E += u5 if ecc < 1.0: nu = E_to_nu(E, ecc) else: if ecc == 1.0: # Parabolic nu = D_to_nu(E) else: # Hyperbolic nu = F_to_nu(E, ecc) return nu
def markley_bulk(k, tof, pp, eecc, iinc, rraan, aargp, nnu0): """ Solves the kepler problem by a non iterative method. Relative error is around 1e-18, only limited by machine double-precission errors. Parameters ---------- k : float Standar Gravitational parameter r0 : array Initial position vector wrt attractor center. v0 : array Initial velocity vector. tof : float Time of flight. Returns ------- rf: array Final position vector vf: array Final velocity vector Note ---- The following algorithm was taken from http://dx.doi.org/10.1007/BF00691917. """ # Solve first for eccentricity and mean anomaly # p, ecc, inc, raan, argp, nu = rv2coe(k, r0, v0) # check of all orbits are elipses # if eecc.max() >= 1: # raise ValueError("All eccentricities must be smaller than 1") EE = nu_to_E(nnu0, eecc) MM0 = E_to_M(EE, eecc) aa = pp / (1 - eecc**2) nn = np.sqrt(k / aa**3) MM = MM0 + nn * tof # Range between -pi and pi MM = (MM + np.pi) % (2 * np.pi) - np.pi # MM = MM % (2 * np.pi) # for i, MM_elmt in enumerate(MM): # if MM_elmt > np.pi: # MM[i] = -(2 * np.pi - MM_elmt) # Equation (20) aalpha = (3 * np.pi**2 + 1.6 * (np.pi - np.abs(MM)) / (1 + eecc)) / (np.pi**2 - 6) # Equation (5) dd = 3 * (1 - eecc) + aalpha * eecc # Equation (9) qq = 2 * aalpha * dd * (1 - eecc) - MM**2 # Equation (10) rr = 3 * aalpha * dd * (dd - 1 + eecc) * MM + MM**3 # Equation (14) ww = (np.abs(rr) + np.sqrt(qq**3 + rr**2))**(2 / 3) # Equation (15) EE = (2 * rr * ww / (ww**2 + ww * qq + qq**2) + MM) / dd # Equation (26) f0 = _kepler_equation(EE, MM, eecc) f1 = _kepler_equation_prime(EE, MM, eecc) f2 = eecc * np.sin(EE) f3 = eecc * np.cos(EE) f4 = -f2 # Equation (22) delta3 = -f0 / (f1 - 0.5 * f0 * f2 / f1) delta4 = -f0 / (f1 + 0.5 * delta3 * f2 + 1 / 6 * delta3**2 * f3) delta5 = -f0 / (f1 + 0.5 * delta4 * f2 + 1 / 6 * delta4**2 * f3 + 1 / 24 * delta4**3 * f4) EE += delta5 nnu = E_to_nu(EE, eecc) return nnu
def rv2coe(k, r, v, tol=1e-8): r"""Converts from vectors to classical orbital elements. 1. First the angular momentum is computed: .. math:: \vec{h} = \vec{r} \times \vec{v} 2. With it the eccentricity can be solved: .. math:: \begin{align} \vec{e} &= \frac{1}{\mu}\left [ \left ( v^{2} - \frac{\mu}{r}\right ) \vec{r} - (\vec{r} \cdot \vec{v})\vec{v} \right ] \\ e &= \sqrt{\vec{e}\cdot\vec{e}} \\ \end{align} 3. The node vector line is solved: .. math:: \begin{align} \vec{N} &= \vec{k} \times \vec{h} \\ N &= \sqrt{\vec{N}\cdot\vec{N}} \end{align} 4. The rigth ascension node is computed: .. math:: \Omega = \left\{ \begin{array}{lcc} cos^{-1}{\left ( \frac{N_{x}}{N} \right )} & if & N_{y} \geq 0 \\ \\ 360^{o} -cos^{-1}{\left ( \frac{N_{x}}{N} \right )} & if & N_{y} < 0 \\ \end{array} \right. 5. The argument of perigee: .. math:: \omega = \left\{ \begin{array}{lcc} cos^{-1}{\left ( \frac{\vec{N}\vec{e}}{Ne} \right )} & if & e_{z} \geq 0 \\ \\ 360^{o} -cos^{-1}{\left ( \frac{\vec{N}\vec{e}}{Ne} \right )} & if & e_{z} < 0 \\ \end{array} \right. 6. And finally the true anomaly: .. math:: \nu = \left\{ \begin{array}{lcc} cos^{-1}{\left ( \frac{\vec{e}\vec{r}}{er} \right )} & if & v_{r} \geq 0 \\ \\ 360^{o} -cos^{-1}{\left ( \frac{\vec{e}\vec{r}}{er} \right )} & if & v_{r} < 0 \\ \end{array} \right. Parameters ---------- k : float Standard gravitational parameter (km^3 / s^2) r : array Position vector (km) v : array Velocity vector (km / s) tol : float, optional Tolerance for eccentricity and inclination checks, default to 1e-8 Returns ------- p : float Semi-latus rectum of parameter (km) ecc: float Eccentricity inc: float Inclination (rad) raan: float Right ascension of the ascending nod (rad) argp: float Argument of Perigee (rad) nu: float True Anomaly (rad) Examples -------- >>> from poliastro.constants import GM_earth >>> from astropy import units as u >>> k = GM_earth.to(u.km ** 3 / u.s ** 2).value # Earth gravitational parameter >>> r = np.array([-6045., -3490., 2500.]) >>> v = np.array([-3.457, 6.618, 2.533]) >>> p, ecc, inc, raan, argp, nu = rv2coe(k, r, v) >>> print("p:", p, "[km]") p: 8530.47436396927 [km] >>> print("ecc:", ecc) ecc: 0.17121118195416898 >>> print("inc:", np.rad2deg(inc), "[deg]") inc: 153.2492285182475 [deg] >>> print("raan:", np.rad2deg(raan), "[deg]") raan: 255.27928533439618 [deg] >>> print("argp:", np.rad2deg(argp), "[deg]") argp: 20.068139973005366 [deg] >>> print("nu:", np.rad2deg(nu), "[deg]") nu: 28.445804984192122 [deg] Note ---- This example is a real exercise from Orbital Mechanics for Engineering students by Howard D.Curtis. This exercise is 4.3 of 3rd. Edition, page 200. """ h = cross(r, v) n = cross([0, 0, 1], h) e = ((v.dot(v) - k / (norm(r))) * r - r.dot(v) * v) / k ecc = norm(e) p = h.dot(h) / k inc = np.arccos(h[2] / norm(h)) circular = ecc < tol equatorial = abs(inc) < tol if equatorial and not circular: raan = 0 argp = np.arctan2(e[1], e[0]) % (2 * np.pi) # Longitude of periapsis nu = np.arctan2(h.dot(cross(e, r)) / norm(h), r.dot(e)) elif not equatorial and circular: raan = np.arctan2(n[1], n[0]) % (2 * np.pi) argp = 0 # Argument of latitude nu = np.arctan2(r.dot(cross(h, n)) / norm(h), r.dot(n)) elif equatorial and circular: raan = 0 argp = 0 nu = np.arctan2(r[1], r[0]) % (2 * np.pi) # True longitude else: a = p / (1 - (ecc**2)) ka = k * a if a > 0: e_se = r.dot(v) / sqrt(ka) e_ce = norm(r) * v.dot(v) / k - 1 nu = E_to_nu(np.arctan2(e_se, e_ce), ecc) else: e_sh = r.dot(v) / sqrt(-ka) e_ch = norm(r) * (norm(v)**2) / k - 1 nu = F_to_nu(np.log((e_ch + e_sh) / (e_ch - e_sh)) / 2, ecc) raan = np.arctan2(n[1], n[0]) % (2 * np.pi) px = r.dot(n) py = r.dot(cross(h, n)) / norm(h) argp = (np.arctan2(py, px) - nu) % (2 * np.pi) nu = (nu + np.pi) % (2 * np.pi) - np.pi return p, ecc, inc, raan, argp, nu
def nu_from_delta_t(delta_t, ecc, k=1.0, q=1.0, delta=1e-2): """True anomaly for given elapsed time since periapsis. Parameters ---------- delta_t : float Time elapsed since periapsis. ecc : float Eccentricity. k : float Gravitational parameter. q : float Periapsis distance. delta : float Parameter that controls the size of the near parabolic region. Returns ------- nu : float True anomaly. """ if ecc < 1 - delta: # Strong elliptic n = np.sqrt(k * (1 - ecc)**3 / q**3) M = n * delta_t # This might represent several revolutions, # so we wrap the true anomaly E = M_to_E((M + np.pi) % (2 * np.pi) - np.pi, ecc) nu = E_to_nu(E, ecc) elif 1 - delta <= ecc < 1: E_delta = np.arccos((1 - delta) / ecc) # We compute M assuming we are in the strong elliptic case # and verify later n = np.sqrt(k * (1 - ecc)**3 / q**3) M = n * delta_t # We check against abs(M) because E_delta could also be negative if E_to_M(E_delta, ecc) <= abs(M): # Strong elliptic, proceed # This might represent several revolutions, # so we wrap the true anomaly E = M_to_E((M + np.pi) % (2 * np.pi) - np.pi, ecc) nu = E_to_nu(E, ecc) else: # Near parabolic, recompute M n = np.sqrt(k / (2 * q**3)) M = n * delta_t D = M_to_D_near_parabolic(M, ecc) nu = D_to_nu(D) elif ecc == 1: # Parabolic n = np.sqrt(k / (2 * q**3)) M = n * delta_t D = M_to_D(M) nu = D_to_nu(D) elif 1 < ecc <= 1 + delta: F_delta = np.arccosh((1 + delta) / ecc) # We compute M assuming we are in the strong hyperbolic case # and verify later n = np.sqrt(k * (ecc - 1)**3 / q**3) M = n * delta_t # We check against abs(M) because F_delta could also be negative if F_to_M(F_delta, ecc) <= abs(M): # Strong hyperbolic, proceed F = M_to_F(M, ecc) nu = F_to_nu(F, ecc) else: # Near parabolic, recompute M n = np.sqrt(k / (2 * q**3)) M = n * delta_t D = M_to_D_near_parabolic(M, ecc) nu = D_to_nu(D) # elif 1 + delta < ecc: else: # Strong hyperbolic n = np.sqrt(k * (ecc - 1)**3 / q**3) M = n * delta_t F = M_to_F(M, ecc) nu = F_to_nu(F, ecc) return nu
def gooding(k, r0, v0, tof, numiter=150, rtol=1e-8): """ Solves the Elliptic Kepler Equation with a cubic convergence and accuracy better than 10e-12 rad is normally achieved. It is not valid for eccentricities equal or higher than 1.0. Parameters ---------- k : float Standard gravitational parameter of the attractor. r : 1x3 vector Position vector. v : 1x3 vector Velocity vector. tof : float Time of flight. rtol: float Relative error for accuracy of the method. Returns ------- rr : 1x3 vector Propagated position vectors. vv : 1x3 vector Note ---- Original paper for the algorithm: https://doi.org/10.1007/BF01238923 """ # Solve first for eccentricity and mean anomaly p, ecc, inc, raan, argp, nu = rv2coe(k, r0, v0) # TODO: parabolic and hyperbolic not implemented cases if ecc >= 1.0: raise NotImplementedError( "Parabolic/Hyperbolic cases still not implemented in gooding.") M0 = nu_to_M(nu, ecc, delta=0) semi_axis_a = p / (1 - ecc**2) n = np.sqrt(k / np.abs(semi_axis_a)**3) M = M0 + n * tof # Start the computation n = 0 c = ecc * np.cos(M) s = ecc * np.sin(M) psi = s / np.sqrt(1 - 2 * c + ecc**2) f = 1.0 while f**2 >= rtol and n <= numiter: xi = np.cos(psi) eta = np.sin(psi) fd = (1 - c * xi) + s * eta fdd = c * eta + s * xi f = psi - fdd psi = psi - f * fd / (fd**2 - 0.5 * f * fdd) n += 1 E = M + psi nu = E_to_nu(E, ecc) return coe2rv(k, p, ecc, inc, raan, argp, nu)
def pimienta(k, r0, v0, tof): """ Raw algorithm for Adonis' Pimienta and John L. Crassidis 15th order polynomial Kepler solver. Parameters ---------- k : float Standar Gravitational parameter r0 : array Initial position vector wrt attractor center. v0 : array Initial velocity vector. tof : float Time of flight. Returns ------- rf: array Final position vector vf: array Final velocity vector Note ---- This algorithm was drived from the original paper: http://hdl.handle.net/10477/50522 """ # TODO: implement hyperbolic case # Solve first for eccentricity and mean anomaly p, ecc, inc, raan, argp, nu = rv2coe(k, r0, v0) M0 = nu_to_M(nu, ecc, delta=0) semi_axis_a = p / (1 - ecc**2) n = np.sqrt(k / np.abs(semi_axis_a)**3) M = M0 + n * tof # Equation (32a), (32b), (32c) and (32d) c3 = 5 / 2 + 560 * ecc a = 15 * (1 - ecc) / c3 b = -M / c3 y = np.sqrt(b**2 / 4 + a**3 / 27) # Equation (33) x_bar = (-b / 2 + y)**(1 / 3) - (b / 2 + y)**(1 / 3) # Coefficients from equations (34a) and (34b) c15 = 3003 / 14336 + 16384 * ecc c13 = 3465 / 13312 - 61440 * ecc c11 = 945 / 2816 + 92160 * ecc c9 = 175 / 384 - 70400 * ecc c7 = 75 / 112 + 28800 * ecc c5 = 9 / 8 - 6048 * ecc # Precompute x_bar powers, equations (35a) to (35d) x_bar2 = x_bar**2 x_bar3 = x_bar2 * x_bar x_bar4 = x_bar3 * x_bar x_bar5 = x_bar4 * x_bar x_bar6 = x_bar5 * x_bar x_bar7 = x_bar6 * x_bar x_bar8 = x_bar7 * x_bar x_bar9 = x_bar8 * x_bar x_bar10 = x_bar9 * x_bar x_bar11 = x_bar10 * x_bar x_bar12 = x_bar11 * x_bar x_bar13 = x_bar12 * x_bar x_bar14 = x_bar13 * x_bar x_bar15 = x_bar14 * x_bar # Function f and its derivatives are given by all the (36) equation set f = (c15 * x_bar15 + c13 * x_bar13 + c11 * x_bar11 + c9 * x_bar9 + c7 * x_bar7 + c5 * x_bar5 + c3 * x_bar3 + 15 * (1 - ecc) * x_bar - M) f1 = (15 * c15 * x_bar14 + 13 * c13 * x_bar12 + 11 * c11 * x_bar10 + 9 * c9 * x_bar8 + 7 * c7 * x_bar6 + 5 * c5 * x_bar4 + 3 * c3 * x_bar2 + 15 * (1 - ecc)) f2 = (210 * c15 * x_bar13 + 156 * c13 * x_bar11 + 110 * c11 * x_bar9 + 72 * c9 * x_bar7 + 42 * c7 * x_bar5 + 20 * c5 * x_bar3 + 6 * c3 * x_bar) f3 = (2730 * c15 * x_bar12 + 1716 * c13 * x_bar10 + 990 * c11 * x_bar8 + 504 * c9 * x_bar6 + 210 * c7 * x_bar4 + 60 * c5 * x_bar2 + 6 * c3) f4 = (32760 * c15 * x_bar11 + 17160 * c13 * x_bar9 + 7920 * c11 * x_bar7 + 3024 * c9 * x_bar5 + 840 * c7 * x_bar3 + 120 * c5 * x_bar) f5 = (360360 * c15 * x_bar10 + 154440 * c13 * x_bar8 + 55440 * c11 * x_bar6 + 15120 * c9 * x_bar4 + 2520 * c7 * x_bar2 + 120 * c5) f6 = (3603600 * c15 * x_bar9 + 1235520 * c13 * x_bar7 + 332640 * c11 * x_bar5 + 60480 * c9 * x_bar3 + 5040 * c7 * x_bar) f7 = (32432400 * c15 * x_bar8 + 8648640 * c13 * x_bar6 + 1663200 * c11 * x_bar4 + 181440 * c9 * x_bar2 + 5040 * c7) f8 = (259459200 * c15 * x_bar7 + 51891840 * c13 * x_bar5 + 6652800 * c11 * x_bar3 + 362880 * c9 * x_bar) f9 = (1.8162144e9 * c15 * x_bar6 + 259459200 * c13 * x_bar4 + 19958400 * c11 * x_bar2 + 362880 * c9) f10 = (1.08972864e10 * c15 * x_bar5 + 1.0378368e9 * c13 * x_bar3 + 39916800 * c11 * x_bar) f11 = 5.4486432e10 * c15 * x_bar4 + 3.1135104e9 * c13 * x_bar2 + 39916800 * c11 f12 = 2.17945728e11 * c15 * x_bar3 + 6.2270208e9 * c13 * x_bar f13 = 6.53837184 * c15 * x_bar2 + 6.2270208e9 * c13 f14 = 1.307674368e13 * c15 * x_bar f15 = 1.307674368e13 * c15 # Solving g parameters defined by equations (37a), (37b), (37c) and (37d) g1 = 1 / 2 g2 = 1 / 6 g3 = 1 / 24 g4 = 1 / 120 g5 = 1 / 720 g6 = 1 / 5040 g7 = 1 / 40320 g8 = 1 / 362880 g9 = 1 / 3628800 g10 = 1 / 39916800 g11 = 1 / 479001600 g12 = 1 / 6.2270208e9 g13 = 1 / 8.71782912e10 g14 = 1 / 1.307674368e12 # Solving for the u_{i} and h_{i} variables defined by equation (38) u1 = -f / f1 h2 = f1 + g1 * u1 * f2 u2 = -f / h2 h3 = f1 + g1 * u2 * f2 + g2 * u2**2 * f3 u3 = -f / h3 h4 = f1 + g1 * u3 * f2 + g2 * u3**2 * f3 + g3 * u3**3 * f4 u4 = -f / h4 h5 = f1 + g1 * u4 * f2 + g2 * u4**2 * f3 + g3 * u4**3 * f4 + g4 * u4**4 * f5 u5 = -f / h5 h6 = (f1 + g1 * u5 * f2 + g2 * u5**2 * f3 + g3 * u5**3 * f4 + g4 * u5**4 * f5 + g5 * u5**5 * f6) u6 = -f / h6 h7 = (f1 + g1 * u6 * f2 + g2 * u6**2 * f3 + g3 * u6**3 * f4 + g4 * u6**4 * f5 + g5 * u6**5 * f6 + g6 * u6**6 * f7) u7 = -f / h7 h8 = (f1 + g1 * u7 * f2 + g2 * u7**2 * f3 + g3 * u7**3 * f4 + g4 * u7**4 * f5 + g5 * u7**5 * f6 + g6 * u7**6 * f7 + g7 * u7**7 * f8) u8 = -f / h8 h9 = (f1 + g1 * u8 * f2 + g2 * u8**2 * f3 + g3 * u8**3 * f4 + g4 * u8**4 * f5 + g5 * u8**5 * f6 + g6 * u8**6 * f7 + g7 * u8**7 * f8 + g8 * u8**8 * f9) u9 = -f / h9 h10 = (f1 + g1 * u9 * f2 + g2 * u9**2 * f3 + g3 * u9**3 * f4 + g4 * u9**4 * f5 + g5 * u9**5 * f6 + g6 * u9**6 * f7 + g7 * u9**7 * f8 + g8 * u9**8 * f9 + g9 * u9**9 * f10) u10 = -f / h10 h11 = (f1 + g1 * u10 * f2 + g2 * u10**2 * f3 + g3 * u10**3 * f4 + g4 * u10**4 * f5 + g5 * u10**5 * f6 + g6 * u10**6 * f7 + g7 * u10**7 * f8 + g8 * u10**8 * f9 + g9 * u10**9 * f10 + g10 * u10**10 * f11) u11 = -f / h11 h12 = (f1 + g1 * u11 * f2 + g2 * u11**2 * f3 + g3 * u11**3 * f4 + g4 * u11**4 * f5 + g5 * u11**5 * f6 + g6 * u11**6 * f7 + g7 * u11**7 * f8 + g8 * u11**8 * f9 + g9 * u11**9 * f10 + g10 * u11**10 * f11 + g11 * u11**11 * f12) u12 = -f / h12 h13 = (f1 + g1 * u12 * f2 + g2 * u12**2 * f3 + g3 * u12**3 * f4 + g4 * u12**4 * f5 + g5 * u12**5 * f6 + g6 * u12**6 * f7 + g7 * u12**7 * f8 + g8 * u12**8 * f9 + g9 * u12**9 * f10 + g10 * u12**10 * f11 + g11 * u12**11 * f12 + g12 * u12**12 * f13) u13 = -f / h13 h14 = (f1 + g1 * u13 * f2 + g2 * u13**2 * f3 + g3 * u13**3 * f4 + g4 * u13**4 * f5 + g5 * u13**5 * f6 + g6 * u13**6 * f7 + g7 * u13**7 * f8 + g8 * u13**8 * f9 + g9 * u13**9 * f10 + g10 * u13**10 * f11 + g11 * u13**11 * f12 + g12 * u13**12 * f13 + g13 * u13**13 * f14) u14 = -f / h14 h15 = (f1 + g1 * u14 * f2 + g2 * u14**2 * f3 + g3 * u14**3 * f4 + g4 * u14**4 * f5 + g5 * u14**5 * f6 + g6 * u14**6 * f7 + g7 * u14**7 * f8 + g8 * u14**8 * f9 + g9 * u14**9 * f10 + g10 * u14**10 * f11 + g11 * u14**11 * f12 + g12 * u14**12 * f13 + g13 * u14**13 * f14 + g14 * u14**14 * f15) u15 = -f / h15 # Solving for x x = x_bar + u15 w = x - 0.01171875 * x**17 / (1 + ecc) # Solving for the true anomaly from eccentricity anomaly E = M + ecc * (-16384 * w**15 + 61440 * w**13 - 92160 * w**11 + 70400 * w**9 - 28800 * w**7 + 6048 * w**5 - 560 * w**3 + 15 * w) nu = E_to_nu(E, ecc) return coe2rv(k, p, ecc, inc, raan, argp, nu)
def markley(k, r0, v0, tof): """ Solves the kepler problem by a non iterative method. Relative error is around 1e-18, only limited by machine double-precission errors. Parameters ---------- k : float Standar Gravitational parameter r0 : array Initial position vector wrt attractor center. v0 : array Initial velocity vector. tof : float Time of flight. Returns ------- rf: array Final position vector vf: array Final velocity vector Note ---- The following algorithm was taken from http://dx.doi.org/10.1007/BF00691917. """ # Solve first for eccentricity and mean anomaly p, ecc, inc, raan, argp, nu = rv2coe(k, r0, v0) M0 = nu_to_M(nu, ecc, delta=0) a = p / (1 - ecc**2) n = np.sqrt(k / a**3) M = M0 + n * tof # Range between -pi and pi M = M % (2 * np.pi) if M > np.pi: M = -(2 * np.pi - M) # Equation (20) alpha = (3 * np.pi**2 + 1.6 * (np.pi - np.abs(M)) / (1 + ecc)) / (np.pi**2 - 6) # Equation (5) d = 3 * (1 - ecc) + alpha * ecc # Equation (9) q = 2 * alpha * d * (1 - ecc) - M**2 # Equation (10) r = 3 * alpha * d * (d - 1 + ecc) * M + M**3 # Equation (14) w = (np.abs(r) + np.sqrt(q**3 + r**2))**(2 / 3) # Equation (15) E = (2 * r * w / (w**2 + w * q + q**2) + M) / d # Equation (26) f0 = _kepler_equation(E, M, ecc) f1 = _kepler_equation_prime(E, M, ecc) f2 = ecc * np.sin(E) f3 = ecc * np.cos(E) f4 = -f2 # Equation (22) delta3 = -f0 / (f1 - 0.5 * f0 * f2 / f1) delta4 = -f0 / (f1 + 0.5 * delta3 * f2 + 1 / 6 * delta3**2 * f3) delta5 = -f0 / (f1 + 0.5 * delta4 * f2 + 1 / 6 * delta4**2 * f3 + 1 / 24 * delta4**3 * f4) E += delta5 nu = E_to_nu(E, ecc) return coe2rv(k, p, ecc, inc, raan, argp, nu)
def mikkola(k, r0, v0, tof, rtol=None): """ Raw algorithm for Mikkola's Kepler solver. Parameters ---------- k : ~astropy.units.Quantity Standard gravitational parameter of the attractor. r : ~astropy.units.Quantity Position vector. v : ~astropy.units.Quantity Velocity vector. tofs : ~astropy.units.Quantity Array of times to propagate. rtol: float This method does not require of tolerance since it is non iterative. Returns ------- rr : ~astropy.units.Quantity Propagated position vectors. vv : ~astropy.units.Quantity Note ---- Original paper: https://doi.org/10.1007/BF01235850 """ # Solving for the classical elements p, ecc, inc, raan, argp, nu = rv2coe(k, r0, v0) M0 = nu_to_M(nu, ecc, delta=0) a = p / (1 - ecc**2) n = np.sqrt(k / np.abs(a)**3) M = M0 + n * tof # Solve for specific geometrical case if ecc < 1.0: # Equation (9a) alpha = (1 - ecc) / (4 * ecc + 1 / 2) else: alpha = (ecc - 1) / (4 * ecc + 1 / 2) beta = M / 2 / (4 * ecc + 1 / 2) # Equation (9b) if beta >= 0: z = (beta + np.sqrt(beta**2 + alpha**3))**(1 / 3) else: z = (beta - np.sqrt(beta**2 + alpha**3))**(1 / 3) s = z - alpha / z # Apply initial correction if ecc < 1.0: ds = -0.078 * s**5 / (1 + ecc) else: ds = 0.071 * s**5 / (1 + 0.45 * s**2) / (1 + 4 * s**2) / ecc s += ds # Solving for the true anomaly if ecc < 1.0: E = M + ecc * (3 * s - 4 * s**3) f = E - ecc * np.sin(E) - M f1 = 1.0 - ecc * np.cos(E) f2 = ecc * np.sin(E) f3 = ecc * np.cos(E) f4 = -f2 f5 = -f3 else: E = 3 * np.log(s + np.sqrt(1 + s**2)) f = -E + ecc * np.sinh(E) - M f1 = -1.0 + ecc * np.cosh(E) f2 = ecc * np.sinh(E) f3 = ecc * np.cosh(E) f4 = f2 f5 = f3 # Apply Taylor expansion u1 = -f / f1 u2 = -f / (f1 + 0.5 * f2 * u1) u3 = -f / (f1 + 0.5 * f2 * u2 + (1.0 / 6.0) * f3 * u2**2) u4 = -f / (f1 + 0.5 * f2 * u3 + (1.0 / 6.0) * f3 * u3**2 + (1.0 / 24.0) * f4 * (u3**3)) u5 = -f / (f1 + f2 * u4 / 2 + f3 * (u4 * u4) / 6.0 + f4 * (u4 * u4 * u4) / 24.0 + f5 * (u4 * u4 * u4 * u4) / 120.0) E += u5 if ecc < 1.0: nu = E_to_nu(E, ecc) else: if ecc == 1.0: # Parabolic nu = D_to_nu(E) else: # Hyperbolic nu = F_to_nu(E, ecc) return coe2rv(k, p, ecc, inc, raan, argp, nu)