Add Universal Variable Lambert solver for computing transfer orbits between any two planets. Enables pork chop plot generation as SQL: SELECT dep_date, arr_date, lambert_c3(3, 4, dep_date, arr_date) FROM generate_series(...) dep CROSS JOIN generate_series(...) arr; New functions: - lambert_transfer(dep_body, arr_body, dep_time, arr_time) → RECORD Returns C3 departure/arrival (km^2/s^2), v_infinity (km/s), time of flight (days), and transfer orbit SMA (AU). - lambert_c3(dep_body, arr_body, dep_time, arr_time) → float8 Convenience: departure C3 only, NULL on solver failure. The solver uses Stumpff functions for unified elliptic/parabolic/hyperbolic handling, with Newton-Raphson iteration and bisection fallback. Each solve is sub-millisecond; PARALLEL SAFE for batch computation. All 11 regression tests pass.
61 lines
1.8 KiB
C
61 lines
1.8 KiB
C
/*
|
|
* lambert.h -- Lambert transfer orbit solver
|
|
*
|
|
* Solves Lambert's problem: given two position vectors and a
|
|
* time of flight, find the orbit connecting them.
|
|
*
|
|
* Uses the Universal Variable formulation with Stumpff functions,
|
|
* handling elliptic, parabolic, and hyperbolic transfers uniformly.
|
|
*
|
|
* Reference:
|
|
* Bate, Mueller & White, "Fundamentals of Astrodynamics" (1971)
|
|
* Battin, "An Introduction to the Methods of Astrodynamics" (1999)
|
|
* Curtis, "Orbital Mechanics for Engineering Students" (2014)
|
|
*
|
|
* All units: AU, days, AU/day. Gravitational parameter is Gauss's
|
|
* constant squared (k^2 = GM_sun in AU^3/day^2).
|
|
*
|
|
* Thread-safe: no static mutable state.
|
|
*/
|
|
|
|
#ifndef PG_ORBIT_LAMBERT_H
|
|
#define PG_ORBIT_LAMBERT_H
|
|
|
|
#ifdef __cplusplus
|
|
extern "C" {
|
|
#endif
|
|
|
|
/*
|
|
* Result of a Lambert transfer orbit solution.
|
|
* All velocities in AU/day.
|
|
*/
|
|
typedef struct lambert_result
|
|
{
|
|
double v1[3]; /* departure velocity vector (AU/day) */
|
|
double v2[3]; /* arrival velocity vector (AU/day) */
|
|
double sma; /* semi-major axis (AU), negative for hyperbolic */
|
|
int converged; /* 1 if solver converged, 0 otherwise */
|
|
} lambert_result;
|
|
|
|
/*
|
|
* Solve Lambert's problem.
|
|
*
|
|
* r1[3]: departure position (AU, heliocentric ecliptic J2000)
|
|
* r2[3]: arrival position (AU, heliocentric ecliptic J2000)
|
|
* tof_days: time of flight in days (must be > 0)
|
|
* mu: gravitational parameter (AU^3/day^2), use GAUSS_K2 for Sun
|
|
* prograde: 1 for prograde (short way), 0 for retrograde (long way)
|
|
* result: output lambert_result
|
|
*
|
|
* Returns 1 on success, 0 on failure.
|
|
*/
|
|
int lambert_solve_uv(const double r1[3], const double r2[3],
|
|
double tof_days, double mu, int prograde,
|
|
lambert_result *result);
|
|
|
|
#ifdef __cplusplus
|
|
}
|
|
#endif
|
|
|
|
#endif /* PG_ORBIT_LAMBERT_H */
|