|
| 1 | +""" A module hosting all algorithms devised by Izzo """ |
| 2 | + |
| 3 | +import numpy as np |
| 4 | +from numpy import cross, pi |
| 5 | +from numpy.linalg import norm |
| 6 | +from scipy.special import hyp2f1 |
| 7 | + |
| 8 | +from lamberthub.utils.assertions import assert_parameters_are_valid |
| 9 | + |
| 10 | + |
| 11 | +def izzo2015( |
| 12 | + mu, |
| 13 | + r1, |
| 14 | + r2, |
| 15 | + tof, |
| 16 | + M=0, |
| 17 | + prograde=True, |
| 18 | + low_path=True, |
| 19 | + maxiter=35, |
| 20 | + atol=1e-5, |
| 21 | + rtol=1e-7, |
| 22 | + full_output=False, |
| 23 | +): |
| 24 | + r""" |
| 25 | + Solves Lambert problem using Izzo's devised algorithm. |
| 26 | +
|
| 27 | + Parameters |
| 28 | + ---------- |
| 29 | + mu: float |
| 30 | + Gravitational parameter, equivalent to :math:`GM` of attractor body. |
| 31 | + r1: numpy.array |
| 32 | + Initial position vector. |
| 33 | + r2: numpy.array |
| 34 | + Final position vector. |
| 35 | + M: int |
| 36 | + Number of revolutions. Must be equal or greater than 0 value. |
| 37 | + prograde: bool |
| 38 | + If `True`, specifies prograde motion. Otherwise, retrograde motion is imposed. |
| 39 | + low_path: bool |
| 40 | + If two solutions are available, it selects between high or low path. |
| 41 | + maxiter: int |
| 42 | + Maximum number of iterations. |
| 43 | + atol: float |
| 44 | + Absolute tolerance. |
| 45 | + rtol: float |
| 46 | + Relative tolerance. |
| 47 | + full_output: bool |
| 48 | + If True, the number of iterations is also returned. |
| 49 | +
|
| 50 | + Returns |
| 51 | + ------- |
| 52 | + v1: numpy.array |
| 53 | + Initial velocity vector. |
| 54 | + v2: numpy.array |
| 55 | + Final velocity vector. |
| 56 | + numiter: list |
| 57 | + Number of iterations. |
| 58 | +
|
| 59 | + Notes |
| 60 | + ----- |
| 61 | + This is the algorithm devised by Dario Izzo[1] in 2015. It inherits from |
| 62 | + the one developed by Lancaster[2] during the 60s, following the universal |
| 63 | + formulae approach. It is one of the most modern solvers, being a complete |
| 64 | + Lambert's problem solver (zero and Multiple-revolution solutions). It shows |
| 65 | + high performance and robustness while requiring no more than four iterations |
| 66 | + to reach a solution. |
| 67 | +
|
| 68 | + All credits of the implementation go to Juan Luis Cano Rodríguez and the |
| 69 | + poliastro development team, from which this routine inherits. Some changes |
| 70 | + were made to adapt it to `lamberthub` API. In addition, the hypergeometric |
| 71 | + function is the one from SciPy. |
| 72 | +
|
| 73 | + Copyright (c) 2012-2021 Juan Luis Cano Rodríguez and the poliastro development team |
| 74 | +
|
| 75 | + References |
| 76 | + ---------- |
| 77 | + [1] Izzo, D. (2015). Revisiting Lambert’s problem. Celestial Mechanics |
| 78 | + and Dynamical Astronomy, 121(1), 1-15. |
| 79 | +
|
| 80 | + [2] Lancaster, E. R., & Blanchard, R. C. (1969). A unified form of |
| 81 | + Lambert's theorem (Vol. 5368). National Aeronautics and Space |
| 82 | + Administration. |
| 83 | +
|
| 84 | + """ |
| 85 | + |
| 86 | + # Check that input parameters are safe |
| 87 | + assert_parameters_are_valid(mu, r1, r2, tof, M) |
| 88 | + |
| 89 | + # Chord |
| 90 | + c = r2 - r1 |
| 91 | + c_norm, r1_norm, r2_norm = norm(c), norm(r1), norm(r2) |
| 92 | + |
| 93 | + # Semiperimeter |
| 94 | + s = (r1_norm + r2_norm + c_norm) * 0.5 |
| 95 | + |
| 96 | + # Versors |
| 97 | + i_r1, i_r2 = r1 / r1_norm, r2 / r2_norm |
| 98 | + i_h = cross(i_r1, i_r2) |
| 99 | + i_h = i_h / norm(i_h) |
| 100 | + |
| 101 | + # Geometry of the problem |
| 102 | + ll = np.sqrt(1 - min(1.0, c_norm / s)) |
| 103 | + |
| 104 | + # Compute the fundamental tangential directions |
| 105 | + if i_h[2] < 0: |
| 106 | + ll = -ll |
| 107 | + i_t1, i_t2 = cross(i_r1, i_h), cross(i_r2, i_h) |
| 108 | + else: |
| 109 | + i_t1, i_t2 = cross(i_h, i_r1), cross(i_h, i_r2) |
| 110 | + |
| 111 | + # Correct transfer angle parameter and tangential vectors regarding orbit's |
| 112 | + # inclination |
| 113 | + ll, i_t1, i_t2 = (-ll, -i_t1, -i_t2) if prograde is False else (ll, i_t1, i_t2) |
| 114 | + |
| 115 | + # Non dimensional time of flight |
| 116 | + T = np.sqrt(2 * mu / s ** 3) * tof |
| 117 | + |
| 118 | + # Find solutions and filter them |
| 119 | + x, y, numiter = _find_xy(ll, T, M, maxiter, rtol, low_path) |
| 120 | + |
| 121 | + # Reconstruct |
| 122 | + gamma = np.sqrt(mu * s / 2) |
| 123 | + rho = (r1_norm - r2_norm) / c_norm |
| 124 | + sigma = np.sqrt(1 - rho ** 2) |
| 125 | + |
| 126 | + # Compute the radial and tangential components at initial and final |
| 127 | + # position vectors |
| 128 | + V_r1, V_r2, V_t1, V_t2 = _reconstruct(x, y, r1_norm, r2_norm, ll, gamma, rho, sigma) |
| 129 | + |
| 130 | + # Solve for the initial and final velocity |
| 131 | + v1 = V_r1 * (r1 / r1_norm) + V_t1 * i_t1 |
| 132 | + v2 = V_r2 * (r2 / r2_norm) + V_t2 * i_t2 |
| 133 | + |
| 134 | + return (v1, v2, numiter) if full_output is True else (v1, v2) |
| 135 | + |
| 136 | + |
| 137 | +def _reconstruct(x, y, r1, r2, ll, gamma, rho, sigma): |
| 138 | + """Reconstruct solution velocity vectors.""" |
| 139 | + V_r1 = gamma * ((ll * y - x) - rho * (ll * y + x)) / r1 |
| 140 | + V_r2 = -gamma * ((ll * y - x) + rho * (ll * y + x)) / r2 |
| 141 | + V_t1 = gamma * sigma * (y + ll * x) / r1 |
| 142 | + V_t2 = gamma * sigma * (y + ll * x) / r2 |
| 143 | + return [V_r1, V_r2, V_t1, V_t2] |
| 144 | + |
| 145 | + |
| 146 | +def _find_xy(ll, T, M, maxiter, rtol, low_path): |
| 147 | + """Computes all x, y for given number of revolutions.""" |
| 148 | + # For abs(ll) == 1 the derivative is not continuous |
| 149 | + assert abs(ll) < 1 |
| 150 | + |
| 151 | + M_max = np.floor(T / pi) |
| 152 | + T_00 = np.arccos(ll) + ll * np.sqrt(1 - ll ** 2) # T_xM |
| 153 | + |
| 154 | + # Refine maximum number of revolutions if necessary |
| 155 | + if T < T_00 + M_max * pi and M_max > 0: |
| 156 | + _, T_min = _compute_T_min(ll, M_max, maxiter, rtol) |
| 157 | + if T < T_min: |
| 158 | + M_max -= 1 |
| 159 | + |
| 160 | + # Check if a feasible solution exist for the given number of revolutions |
| 161 | + # This departs from the original paper in that we do not compute all solutions |
| 162 | + if M > M_max: |
| 163 | + raise ValueError("No feasible solution, try lower M!") |
| 164 | + |
| 165 | + # Initial guess |
| 166 | + x_0 = _initial_guess(T, ll, M, low_path) |
| 167 | + |
| 168 | + # Start Householder iterations from x_0 and find x, y |
| 169 | + x, numiter = _householder(x_0, T, ll, M, rtol, maxiter) |
| 170 | + y = _compute_y(x, ll) |
| 171 | + |
| 172 | + return x, y, numiter |
| 173 | + |
| 174 | + |
| 175 | +def _compute_y(x, ll): |
| 176 | + """Computes y.""" |
| 177 | + return np.sqrt(1 - ll ** 2 * (1 - x ** 2)) |
| 178 | + |
| 179 | + |
| 180 | +def _compute_psi(x, y, ll): |
| 181 | + """Computes psi. |
| 182 | +
|
| 183 | + "The auxiliary angle psi is computed using Eq.(17) by the appropriate |
| 184 | + inverse function" |
| 185 | +
|
| 186 | + """ |
| 187 | + if -1 <= x < 1: |
| 188 | + # Elliptic motion |
| 189 | + # Use arc cosine to avoid numerical errors |
| 190 | + return np.arccos(x * y + ll * (1 - x ** 2)) |
| 191 | + elif x > 1: |
| 192 | + # Hyperbolic motion |
| 193 | + # The hyperbolic sine is bijective |
| 194 | + return np.arcsinh((y - x * ll) * np.sqrt(x ** 2 - 1)) |
| 195 | + else: |
| 196 | + # Parabolic motion |
| 197 | + return 0.0 |
| 198 | + |
| 199 | + |
| 200 | +def _tof_equation(x, T0, ll, M): |
| 201 | + """Time of flight equation.""" |
| 202 | + return _tof_equation_y(x, _compute_y(x, ll), T0, ll, M) |
| 203 | + |
| 204 | + |
| 205 | +def _tof_equation_y(x, y, T0, ll, M): |
| 206 | + """Time of flight equation with externally computated y.""" |
| 207 | + if M == 0 and np.sqrt(0.6) < x < np.sqrt(1.4): |
| 208 | + eta = y - ll * x |
| 209 | + S_1 = (1 - ll - x * eta) * 0.5 |
| 210 | + Q = 4 / 3 * hyp2f1(3, 1, 5 / 2, S_1) |
| 211 | + T_ = (eta ** 3 * Q + 4 * ll * eta) * 0.5 |
| 212 | + else: |
| 213 | + psi = _compute_psi(x, y, ll) |
| 214 | + T_ = np.divide( |
| 215 | + np.divide(psi + M * pi, np.sqrt(np.abs(1 - x ** 2))) - x + ll * y, |
| 216 | + (1 - x ** 2), |
| 217 | + ) |
| 218 | + |
| 219 | + return T_ - T0 |
| 220 | + |
| 221 | + |
| 222 | +def _tof_equation_p(x, y, T, ll): |
| 223 | + # TODO: What about derivatives when x approaches 1? |
| 224 | + return (3 * T * x - 2 + 2 * ll ** 3 * x / y) / (1 - x ** 2) |
| 225 | + |
| 226 | + |
| 227 | +def _tof_equation_p2(x, y, T, dT, ll): |
| 228 | + return (3 * T + 5 * x * dT + 2 * (1 - ll ** 2) * ll ** 3 / y ** 3) / (1 - x ** 2) |
| 229 | + |
| 230 | + |
| 231 | +def _tof_equation_p3(x, y, _, dT, ddT, ll): |
| 232 | + return (7 * x * ddT + 8 * dT - 6 * (1 - ll ** 2) * ll ** 5 * x / y ** 5) / ( |
| 233 | + 1 - x ** 2 |
| 234 | + ) |
| 235 | + |
| 236 | + |
| 237 | +def _compute_T_min(ll, M, maxiter, rtol): |
| 238 | + """Compute minimum T.""" |
| 239 | + if ll == 1: |
| 240 | + x_T_min = 0.0 |
| 241 | + T_min = _tof_equation(x_T_min, 0.0, ll, M) |
| 242 | + else: |
| 243 | + if M == 0: |
| 244 | + x_T_min = np.inf |
| 245 | + T_min = 0.0 |
| 246 | + else: |
| 247 | + # Set x_i > 0 to avoid problems at ll = -1 |
| 248 | + x_i = 0.1 |
| 249 | + T_i = _tof_equation(x_i, 0.0, ll, M) |
| 250 | + x_T_min = _halley(x_i, T_i, ll, rtol, maxiter) |
| 251 | + T_min = _tof_equation(x_T_min, 0.0, ll, M) |
| 252 | + |
| 253 | + return [x_T_min, T_min] |
| 254 | + |
| 255 | + |
| 256 | +def _initial_guess(T, ll, M, low_path): |
| 257 | + """Initial guess.""" |
| 258 | + if M == 0: |
| 259 | + # Single revolution |
| 260 | + T_0 = np.arccos(ll) + ll * np.sqrt(1 - ll ** 2) + M * pi # Equation 19 |
| 261 | + T_1 = 2 * (1 - ll ** 3) / 3 # Equation 21 |
| 262 | + if T >= T_0: |
| 263 | + x_0 = (T_0 / T) ** (2 / 3) - 1 |
| 264 | + elif T < T_1: |
| 265 | + x_0 = 5 / 2 * T_1 / T * (T_1 - T) / (1 - ll ** 5) + 1 |
| 266 | + else: |
| 267 | + # This is the real condition, which is not exactly equivalent |
| 268 | + # elif T_1 < T < T_0 |
| 269 | + x_0 = (T_0 / T) ** (np.log2(T_1 / T_0)) - 1 |
| 270 | + |
| 271 | + return x_0 |
| 272 | + else: |
| 273 | + # Multiple revolution |
| 274 | + x_0l = (((M * pi + pi) / (8 * T)) ** (2 / 3) - 1) / ( |
| 275 | + ((M * pi + pi) / (8 * T)) ** (2 / 3) + 1 |
| 276 | + ) |
| 277 | + x_0r = (((8 * T) / (M * pi)) ** (2 / 3) - 1) / ( |
| 278 | + ((8 * T) / (M * pi)) ** (2 / 3) + 1 |
| 279 | + ) |
| 280 | + |
| 281 | + # Filter out the solution |
| 282 | + x_0 = np.max([x_0l, x_0r]) if low_path is True else np.min([x_0l, x_0r]) |
| 283 | + |
| 284 | + return x_0 |
| 285 | + |
| 286 | + |
| 287 | +def _halley(p0, T0, ll, tol, maxiter): |
| 288 | + """Find a minimum of time of flight equation using the Halley method. |
| 289 | +
|
| 290 | + Note |
| 291 | + ---- |
| 292 | + This function is private because it assumes a calling convention specific to |
| 293 | + this module and is not really reusable. |
| 294 | +
|
| 295 | + """ |
| 296 | + for ii in range(1, maxiter + 1): |
| 297 | + y = _compute_y(p0, ll) |
| 298 | + fder = _tof_equation_p(p0, y, T0, ll) |
| 299 | + fder2 = _tof_equation_p2(p0, y, T0, fder, ll) |
| 300 | + if fder2 == 0: |
| 301 | + raise RuntimeError("Derivative was zero") |
| 302 | + fder3 = _tof_equation_p3(p0, y, T0, fder, fder2, ll) |
| 303 | + |
| 304 | + # Halley step (cubic) |
| 305 | + p = p0 - 2 * fder * fder2 / (2 * fder2 ** 2 - fder * fder3) |
| 306 | + |
| 307 | + if abs(p - p0) < tol: |
| 308 | + return p |
| 309 | + p0 = p |
| 310 | + |
| 311 | + raise RuntimeError("Failed to converge") |
| 312 | + |
| 313 | + |
| 314 | +def _householder(p0, T0, ll, M, tol, maxiter): |
| 315 | + """Find a zero of time of flight equation using the Householder method. |
| 316 | +
|
| 317 | + Note |
| 318 | + ---- |
| 319 | + This function is private because it assumes a calling convention specific to |
| 320 | + this module and is not really reusable. |
| 321 | +
|
| 322 | + """ |
| 323 | + for numiter in range(1, maxiter + 1): |
| 324 | + y = _compute_y(p0, ll) |
| 325 | + fval = _tof_equation_y(p0, y, T0, ll, M) |
| 326 | + T = fval + T0 |
| 327 | + fder = _tof_equation_p(p0, y, T, ll) |
| 328 | + fder2 = _tof_equation_p2(p0, y, T, fder, ll) |
| 329 | + fder3 = _tof_equation_p3(p0, y, T, fder, fder2, ll) |
| 330 | + |
| 331 | + # Householder step (quartic) |
| 332 | + p = p0 - fval * ( |
| 333 | + (fder ** 2 - fval * fder2 / 2) |
| 334 | + / (fder * (fder ** 2 - fval * fder2) + fder3 * fval ** 2 / 6) |
| 335 | + ) |
| 336 | + |
| 337 | + if abs(p - p0) < tol: |
| 338 | + return p, numiter |
| 339 | + p0 = p |
| 340 | + |
| 341 | + raise RuntimeError("Failed to converge") |
0 commit comments