pg_orrery/src/od_solver.h
Ryan Malloy adfb6949e1 Add range rate fitting, weighted observations, and Gauss angles-only IOD (v0.6.0)
Range rate: topocentric residuals now include an optional 4th component
(dot(Δr, v_ecef) / |Δr|) with OD_RR_SCALE=10.0 for unit balancing.
Controlled via fit_range_rate parameter on tle_from_topocentric().

Weighted observations: per-observation weights applied as √w scaling
to both residuals and Jacobian rows, producing the weighted normal
equations H'WH without explicit W construction. Weights parameter
added to tle_from_eci, tle_from_topocentric, and tle_from_angles.

Gauss angles-only IOD: Vallado Algorithm 52 implementation for
seed-free orbit recovery from 3+ RA/Dec observations. New RA/Dec
residual function with cos(dec) scaling and wrap-around handling.
New tle_from_angles() and tle_from_angles_multi() SQL functions
accepting RA in hours [0,24), Dec in degrees [-90,90].

New standalone test suite: test_od_gauss (17 assertions).
New regression tests: Tests 18-25 covering range rate, weights,
angles-only with/without seed, and error cases.
2026-02-17 17:48:13 -07:00

110 lines
3.5 KiB
C

/*
* od_solver.h -- TLE fitting via weighted least-squares differential correction
*
* Implements Vallado & Crawford (2008) AIAA 2008-6770:
* - Equinoctial element perturbation
* - SGP4 as the propagation engine
* - Jacobian by finite differencing
* - SVD solve via LAPACK dgelss_()
* - Tiered step limiting (Vallado 2008 Section V.B)
* - Brouwer-to-Kozai Newton-Raphson for TLE output
*/
#ifndef PG_ORRERY_OD_SOLVER_H
#define PG_ORRERY_OD_SOLVER_H
#include "norad.h"
/* Maximum number of DC iterations */
#define OD_MAX_ITER 25
/* Convergence thresholds (Vallado 2008 Section V.A) */
#define OD_RMS_REL_TOL 0.0002
#define OD_RMS_ABS_TOL 0.0002
/* Number of equinoctial states (6 orbital + optional B*) */
#define OD_NSTATE_6 6
#define OD_NSTATE_7 7
/*
* Observation type flag
*/
typedef enum
{
OD_OBS_ECI = 0, /* 6-component: x, y, z, vx, vy, vz (km, km/s) */
OD_OBS_TOPO = 1, /* 3-component: az, el, range (rad, rad, km) */
OD_OBS_RADEC = 2 /* 2-component: RA, Dec (radians, equatorial) */
} od_obs_type_t;
/*
* Single observation for the DC solver
*/
typedef struct
{
double jd; /* Julian date of observation */
double data[6]; /* ECI: [x,y,z,vx,vy,vz], topo: [az,el,range,...] */
int observer_idx; /* index into config->observers[] (topo mode) */
} od_observation_t;
/*
* Observer location (needed for topocentric mode)
*/
typedef struct
{
double lat; /* radians */
double lon; /* radians */
double alt_m; /* meters */
double ecef[3]; /* pre-computed ECEF position (km) */
} od_observer_t;
/*
* Solver configuration
*/
typedef struct
{
od_obs_type_t obs_type; /* ECI or topocentric */
int fit_bstar; /* include B* as 7th state */
int fit_range_rate; /* include range_rate in topo residuals */
int max_iter; /* iteration limit */
od_observer_t *observers; /* array of observers (topo mode) */
int n_observers; /* count (0 for ECI mode) */
double *weights; /* per-observation weights (NULL=uniform) */
int n_weights;
} od_config_t;
/*
* Solver result
*/
typedef struct
{
tle_t fitted_tle; /* resulting TLE (Kozai mean motion) */
int iterations; /* iterations performed */
double rms_initial; /* RMS residual before fitting (km) */
double rms_final; /* RMS residual after fitting (km) */
int converged; /* 1 if converged, 0 if hit max_iter */
char status[64]; /* human-readable status */
double covariance[28]; /* upper triangle of (H^T H)^{-1}, row-major */
int cov_size; /* 21 (6-state) or 28 (7-state), 0 if failed */
double condition_number; /* s_max / s_min from final SVD */
} od_result_t;
/*
* Primary entry point: fit a TLE from observations.
*
* obs: array of N observations
* n_obs: number of observations (must be >= 6, or >= 7 if fit_bstar)
* seed: initial TLE guess (if NULL, computed from first/last obs)
* config: solver configuration
* result: output (caller allocates)
*
* Returns 0 on success, -1 on error.
*
* Memory: all working arrays allocated via palloc() (when used from
* PostgreSQL) or malloc() (standalone). Freed before return.
*/
int od_fit_tle(const od_observation_t *obs, int n_obs,
const tle_t *seed, const od_config_t *config,
od_result_t *result);
#endif /* PG_ORRERY_OD_SOLVER_H */