Skip to main content

ballistics_engine/
trajectory_integration.rs

1//! Advanced trajectory integration methods (RK4, RK45)
2//!
3//! This module provides production-grade numerical integration for ballistic trajectories:
4//! - RK4: 4th-order Runge-Kutta (fixed step)
5//! - RK45: Dormand-Prince adaptive method (same as scipy.integrate.solve_ivp)
6//!
7//! MBA-155: Upstreamed from ballistics_rust for shared use
8
9use nalgebra::{Vector3, Vector6};
10use std::collections::HashMap;
11
12use crate::derivatives::compute_derivatives;
13use crate::wind::WindSegment;
14use crate::BallisticInputs;
15use crate::DragModel;
16
17/// RK4 integration step
18fn rk4_step(
19    state: &Vector6<f64>,
20    t: f64,
21    dt: f64,
22    params: &TrajectoryParams,
23    inputs: &BallisticInputs,
24) -> Vector6<f64> {
25    // RK4 integration
26    let k1 = compute_derivatives_vec(state, t, params, inputs);
27    let k2 = compute_derivatives_vec(&(state + dt * 0.5 * k1), t + dt * 0.5, params, inputs);
28    let k3 = compute_derivatives_vec(&(state + dt * 0.5 * k2), t + dt * 0.5, params, inputs);
29    let k4 = compute_derivatives_vec(&(state + dt * k3), t + dt, params, inputs);
30
31    state + (dt / 6.0) * (k1 + 2.0 * k2 + 2.0 * k3 + k4)
32}
33
34/// Adaptive RK45 integration step (Dormand-Prince method)
35fn rk45_step(
36    state: &Vector6<f64>,
37    t: f64,
38    dt: f64,
39    params: &TrajectoryParams,
40    inputs: &BallisticInputs,
41    tol: f64,
42) -> (Vector6<f64>, f64, f64) {
43    // Dormand-Prince coefficients (same as scipy.integrate.solve_ivp RK45)
44    const A21: f64 = 1.0 / 5.0;
45    const A31: f64 = 3.0 / 40.0;
46    const A32: f64 = 9.0 / 40.0;
47    const A41: f64 = 44.0 / 45.0;
48    const A42: f64 = -56.0 / 15.0;
49    const A43: f64 = 32.0 / 9.0;
50    const A51: f64 = 19372.0 / 6561.0;
51    const A52: f64 = -25360.0 / 2187.0;
52    const A53: f64 = 64448.0 / 6561.0;
53    const A54: f64 = -212.0 / 729.0;
54    const A61: f64 = 9017.0 / 3168.0;
55    const A62: f64 = -355.0 / 33.0;
56    const A63: f64 = 46732.0 / 5247.0;
57    const A64: f64 = 49.0 / 176.0;
58    const A65: f64 = -5103.0 / 18656.0;
59    const A71: f64 = 35.0 / 384.0;
60    const A73: f64 = 500.0 / 1113.0;
61    const A74: f64 = 125.0 / 192.0;
62    const A75: f64 = -2187.0 / 6784.0;
63    const A76: f64 = 11.0 / 84.0;
64
65    // 5th order coefficients
66    const B1: f64 = 35.0 / 384.0;
67    const B3: f64 = 500.0 / 1113.0;
68    const B4: f64 = 125.0 / 192.0;
69    const B5: f64 = -2187.0 / 6784.0;
70    const B6: f64 = 11.0 / 84.0;
71
72    // 4th order coefficients (for error estimation)
73    const B1_ERR: f64 = 5179.0 / 57600.0;
74    const B3_ERR: f64 = 7571.0 / 16695.0;
75    const B4_ERR: f64 = 393.0 / 640.0;
76    const B5_ERR: f64 = -92097.0 / 339200.0;
77    const B6_ERR: f64 = 187.0 / 2100.0;
78    const B7_ERR: f64 = 1.0 / 40.0;
79
80    // Compute stages
81    let k1 = compute_derivatives_vec(state, t, params, inputs);
82    let k2 = compute_derivatives_vec(&(state + dt * A21 * k1), t + dt * 0.2, params, inputs);
83    let k3 =
84        compute_derivatives_vec(&(state + dt * (A31 * k1 + A32 * k2)), t + dt * 0.3, params, inputs);
85    let k4 = compute_derivatives_vec(
86        &(state + dt * (A41 * k1 + A42 * k2 + A43 * k3)),
87        t + dt * 0.8,
88        params,
89        inputs,
90    );
91    let k5 = compute_derivatives_vec(
92        &(state + dt * (A51 * k1 + A52 * k2 + A53 * k3 + A54 * k4)),
93        t + dt * 8.0 / 9.0,
94        params,
95        inputs,
96    );
97    let k6 = compute_derivatives_vec(
98        &(state + dt * (A61 * k1 + A62 * k2 + A63 * k3 + A64 * k4 + A65 * k5)),
99        t + dt,
100        params,
101        inputs,
102    );
103    let k7 = compute_derivatives_vec(
104        &(state + dt * (A71 * k1 + A73 * k3 + A74 * k4 + A75 * k5 + A76 * k6)),
105        t + dt,
106        params,
107        inputs,
108    );
109
110    // 5th order solution
111    let y_new = state + dt * (B1 * k1 + B3 * k3 + B4 * k4 + B5 * k5 + B6 * k6);
112
113    // 4th order solution for error estimate
114    let y_err = state
115        + dt * (B1_ERR * k1 + B3_ERR * k3 + B4_ERR * k4 + B5_ERR * k5 + B6_ERR * k6 + B7_ERR * k7);
116
117    // Error estimate
118    let error = (y_new - y_err).norm() / (1.0 + state.norm());
119
120    // Adaptive step size
121    let safety = 0.9;
122    let dt_new = if error < tol {
123        dt * safety * (tol / error).powf(0.2).min(2.0)
124    } else {
125        dt * safety * (tol / error).powf(0.25).max(0.1)
126    };
127
128    (y_new, dt_new, error)
129}
130
131/// Parameters for trajectory computation
132pub struct TrajectoryParams {
133    pub mass_kg: f64,
134    pub bc: f64,
135    pub drag_model: DragModel,
136    pub wind_segments: Vec<WindSegment>,
137    pub atmos_params: (f64, f64, f64, f64),
138    pub omega_vector: Option<Vector3<f64>>,
139    pub enable_spin_drift: bool,
140    pub enable_magnus: bool,
141    pub enable_coriolis: bool,
142    pub target_distance_m: f64, // Target horizontal distance in meters
143    pub enable_wind_shear: bool,
144    pub wind_shear_model: String,
145    pub shooter_altitude_m: f64,
146    pub is_twist_right: bool, // True for right-hand twist, false for left-hand
147    pub custom_drag_table: Option<crate::drag::DragTable>, // Custom Drag Model (CDM) data
148    pub bc_segments: Option<Vec<(f64, f64)>>, // Mach-based BC segments: (mach, bc)
149    pub use_bc_segments: bool, // Whether to use BC segment interpolation
150    /// MBA-954: altitude (m, relative to launch) below which integration stops. -1000.0 is the
151    /// historical default — effectively "no early ground impact" for normal flat-fire shots.
152    pub ground_threshold: f64,
153}
154
155/// Build the loop-invariant BallisticInputs for the derivatives function ONCE per integration,
156/// instead of rebuilding it (a "none".to_string() alloc plus bc_segments / custom_drag_table
157/// clones) on every derivative evaluation (4x per RK4 step, 7x per RK45 step). muzzle_velocity is
158/// set to 0.0 because compute_derivatives never reads it on this path; every other field depends
159/// only on `params`, so the struct is constant for the whole integration.
160fn build_inputs(params: &TrajectoryParams) -> BallisticInputs {
161    let mut inputs = BallisticInputs {
162        bc_value: params.bc,
163        bc_type: params.drag_model,
164        bullet_mass: params.mass_kg, // kg
165        muzzle_velocity: 0.0,        // unread by compute_derivatives on this path
166        bullet_diameter: 0.0078232,  // 0.308 in -> meters
167        bullet_length: 0.031496,     // 1.24 in -> meters
168        twist_rate: 10.0,            // inches/turn (twist stays imperial)
169        is_twist_right: params.is_twist_right,
170        enable_advanced_effects: params.enable_spin_drift
171            || params.enable_magnus
172            || params.enable_coriolis,
173        enable_magnus: params.enable_magnus,
174        enable_coriolis: params.enable_coriolis,
175        altitude: params.atmos_params.0,
176        temperature: params.atmos_params.1,
177        pressure: params.atmos_params.2,
178        humidity: params.atmos_params.3,
179        tipoff_yaw: 0.0,
180        target_distance: 1000.0, // default
181        muzzle_angle: 0.0,
182        wind_speed: if !params.wind_segments.is_empty() {
183            params.wind_segments[0].0 * 0.2777778 // km/h -> m/s
184        } else {
185            0.0
186        },
187        wind_angle: if !params.wind_segments.is_empty() {
188            params.wind_segments[0].1.to_radians() // degrees -> radians
189        } else {
190            0.0
191        },
192        latitude: None,
193        shooting_angle: 0.0,
194        azimuth_angle: 0.0,
195        use_powder_sensitivity: false,
196        powder_temp_sensitivity: 0.0,
197        powder_temp: 59.0,
198        tipoff_decay_distance: 0.0,
199        ground_threshold: params.ground_threshold, // MBA-954: honor the configured ground plane
200        bc_segments: params.bc_segments.clone(),
201        caliber_inches: 0.308,
202        weight_grains: params.mass_kg / 0.00006479891,
203        use_bc_segments: params.use_bc_segments,
204        bullet_id: None,
205        bc_segments_data: None,
206        use_enhanced_spin_drift: params.enable_spin_drift,
207        use_form_factor: false,
208        manufacturer: None,
209        bullet_model: None,
210        enable_wind_shear: false,
211        wind_shear_model: "none".to_string(),
212        use_cluster_bc: false,
213        bullet_cluster: None,
214        custom_drag_table: params.custom_drag_table.clone(),
215        bc_type_str: None,
216        enable_pitch_damping: false,
217        enable_precession_nutation: false,
218        // MBA-959: aerodynamic jump is intentionally OFF on this path. integrate_trajectory
219        // is a low-level state integrator: it advances a raw initial_state (no muzzle angle to
220        // perturb), carries placeholder bullet geometry (fixed twist/length/diameter) and an
221        // unread muzzle_velocity, so a meaningful Litz Sg can't be formed here (Sg would be ~0
222        // and the AJ guard would suppress it regardless). AJ belongs on the BallisticInputs +
223        // TrajectorySolver path, which bindings use and which already honors the flag.
224        enable_aerodynamic_jump: false,
225        use_rk4: true,
226        use_adaptive_rk45: false,
227        enable_trajectory_sampling: false,
228        sample_interval: 10.0,
229        sight_height: 0.0,
230        muzzle_height: 0.0,
231        target_height: 0.0,
232    };
233
234    // MBA-955: pre-populate velocity-BC segments ONCE here, instead of get_bc_for_velocity
235    // rebuilding them (a model String + a segment Vec) on every derivative evaluation (4-7x per
236    // step). Gated to EXACTLY the case where the per-step path would estimate: use_bc_segments on,
237    // no explicit velocity segments, and no Mach-based bc_segments (those take a different,
238    // unchanged path). bc_used there == params.bc == inputs.bc_value, so the estimated segments are
239    // identical and the per-step fast-path lookup returns the same BC -> byte-identical output.
240    if inputs.use_bc_segments && inputs.bc_segments_data.is_none() && inputs.bc_segments.is_none() {
241        inputs.bc_segments_data =
242            crate::derivatives::estimate_bc_segments_for(&inputs, inputs.bc_value);
243    }
244    inputs
245}
246
247/// Convert state to Vector6 and call compute_derivatives
248fn compute_derivatives_vec(
249    state: &Vector6<f64>,
250    t: f64,
251    params: &TrajectoryParams,
252    inputs: &BallisticInputs,
253) -> Vector6<f64> {
254    let pos = Vector3::new(state[0], state[1], state[2]);
255    let vel = Vector3::new(state[3], state[4], state[5]);
256
257    // Calculate wind at current position with shear support
258    let wind_vector = if !params.wind_segments.is_empty() {
259        if params.enable_wind_shear && params.wind_shear_model != "none" {
260            crate::wind_shear::get_wind_at_position(
261                &pos,
262                &params.wind_segments,
263                params.enable_wind_shear,
264                &params.wind_shear_model,
265                params.shooter_altitude_m,
266            )
267        } else {
268            // Simple constant wind (original implementation)
269            let seg = &params.wind_segments[0];
270            let wind_speed_mps = seg.0 * 0.2777778; // km/h to m/s
271            let wind_angle_rad = seg.1.to_radians();
272            // McCoy: x=downrange, y=vertical, z=lateral
273            Vector3::new(
274                -wind_speed_mps * wind_angle_rad.cos(), // x (downrange - head/tail component)
275                0.0,                                    // y (vertical)
276                -wind_speed_mps * wind_angle_rad.sin(), // z (lateral - crosswind component)
277            )
278        }
279    } else {
280        Vector3::zeros()
281    };
282
283    // Call compute_derivatives - returns [f64; 6] directly. `inputs` is built once per
284    // integration by build_inputs() and threaded in, instead of rebuilt every call.
285    let deriv_result = compute_derivatives(
286        pos,
287        vel,
288        inputs,
289        wind_vector,
290        params.atmos_params,
291        params.bc,
292        params.omega_vector,
293        t,
294    );
295
296    Vector6::new(
297        deriv_result[0],
298        deriv_result[1],
299        deriv_result[2],
300        deriv_result[3],
301        deriv_result[4],
302        deriv_result[5],
303    )
304}
305
306/// Main trajectory integration function
307pub fn integrate_trajectory(
308    initial_state: [f64; 6],
309    t_span: (f64, f64),
310    params: TrajectoryParams,
311    method: &str,
312    tolerance: f64,
313    max_step: f64,
314) -> Vec<(f64, Vector6<f64>)> {
315    let mut state = Vector6::new(
316        initial_state[0],
317        initial_state[1],
318        initial_state[2],
319        initial_state[3],
320        initial_state[4],
321        initial_state[5],
322    );
323
324    let mut t = t_span.0;
325    let t_end = t_span.1;
326    let mut dt = (t_end - t) / 1000.0; // Initial step size
327
328    let mut trajectory = Vec::with_capacity(10000);
329    trajectory.push((t, state));
330
331    // Build the (loop-invariant) derivative inputs once for the whole integration, instead of
332    // rebuilding the struct on every derivative evaluation.
333    let inputs = build_inputs(&params);
334
335    match method {
336        "RK4" => {
337            // Fixed step RK4 with target detection
338            dt = dt.min(max_step).min(0.001); // Use smaller steps for accuracy
339
340            while t < t_end {
341                if t + dt > t_end {
342                    dt = t_end - t;
343                }
344
345                let new_state = rk4_step(&state, t, dt, &params, &inputs);
346
347                // Check if we're about to pass the target (X is downrange, McCoy)
348                if state[0] < params.target_distance_m && new_state[0] >= params.target_distance_m {
349                    // Interpolate to find exact target crossing
350                    let alpha = (params.target_distance_m - state[0]) / (new_state[0] - state[0]);
351                    let dt_to_target = dt * alpha;
352
353                    // Take a smaller step to reach target exactly
354                    let final_state = rk4_step(&state, t, dt_to_target, &params, &inputs);
355
356                    // Ensure we don't overshoot
357                    let mut corrected_state = final_state;
358                    if corrected_state[0] > params.target_distance_m {
359                        corrected_state[0] = params.target_distance_m;
360                    }
361
362                    trajectory.push((t + dt_to_target, corrected_state));
363                    break; // Stop at target
364                }
365
366                state = new_state;
367                t += dt;
368                trajectory.push((t, state));
369
370                // Check if we've reached or passed the target
371                if state[0] >= params.target_distance_m {
372                    // X is downrange (McCoy)
373                    // Add final point exactly at target
374                    let mut final_state = state;
375                    final_state[0] = params.target_distance_m; // X is downrange
376                    trajectory.push((t, final_state));
377                    break;
378                }
379
380                // Check if bullet hit ground (MBA-954: honor the configured ground plane,
381                // not a hardcoded -1000.0)
382                if state[1] < params.ground_threshold {
383                    break;
384                }
385            }
386        }
387        "RK45" | _ => {
388            // Adaptive RK45 with better sampling
389            let mut last_save_x = 0.0; // X is downrange (McCoy)
390            let save_interval_m = params.target_distance_m / 50.0; // Save ~50 points minimum
391
392            // OPTIMIZATION: Adjust max step size when wind shear is enabled
393            // This improves numerical stability at long ranges
394            let effective_max_step =
395                if params.enable_wind_shear && params.wind_shear_model != "none" {
396                    // Use smaller steps for wind shear, but not TOO small
397                    if params.target_distance_m > 800.0 {
398                        0.01 // Smaller steps for long range with shear (10ms)
399                    } else {
400                        0.02 // Normal steps for medium range with shear (20ms)
401                    }
402                } else {
403                    max_step // Use provided max_step when no wind shear
404                };
405
406            // Set initial step size - ensure it's reasonable
407            dt = dt.min(effective_max_step).max(0.0001); // At least 0.1ms to avoid infinite loops
408
409            // Safety check: maximum iterations to prevent infinite loops
410            let max_iterations = 100000; // Should be more than enough for any realistic trajectory
411            let mut iteration_count = 0;
412
413            while t < t_end && iteration_count < max_iterations {
414                iteration_count += 1;
415
416                // Limit time step for better resolution
417                if t + dt > t_end {
418                    dt = t_end - t;
419                }
420
421                let (new_state, dt_new, _error) =
422                    rk45_step(&state, t, dt, &params, &inputs, tolerance);
423
424                // Check if we're about to pass the target (X is downrange, McCoy)
425                if state[0] < params.target_distance_m && new_state[0] >= params.target_distance_m {
426                    // Interpolate to find exact target crossing
427                    let alpha = (params.target_distance_m - state[0]) / (new_state[0] - state[0]);
428                    let dt_to_target = dt * alpha;
429
430                    // Take a smaller step to reach target exactly
431                    let (final_state, _, _) =
432                        rk45_step(&state, t, dt_to_target, &params, &inputs, tolerance);
433
434                    // Make sure we don't overshoot
435                    let mut corrected_state = final_state;
436                    if corrected_state[0] > params.target_distance_m {
437                        corrected_state[0] = params.target_distance_m;
438                    }
439
440                    trajectory.push((t + dt_to_target, corrected_state));
441                    break; // Stop at target - no more points after this
442                }
443
444                // Update state
445                state = new_state;
446                t += dt;
447
448                // Save trajectory point if we've moved enough distance
449                if state[0] - last_save_x >= save_interval_m || state[0] >= params.target_distance_m
450                {
451                    // X is downrange
452                    trajectory.push((t, state));
453                    last_save_x = state[0];
454                }
455
456                // Limit dt for next step - ensure we get enough resolution
457                dt = dt_new.min(effective_max_step).max(0.0001); // Use effective max step, min 0.1ms
458
459                // Stop if we've reached the target
460                if state[0] >= params.target_distance_m {
461                    // X is downrange (McCoy)
462                    // Add final point at target distance
463                    let mut final_state = state;
464                    final_state[0] = params.target_distance_m; // X is downrange
465                    trajectory.push((t, final_state));
466                    break;
467                }
468
469                // Check if bullet hit ground (MBA-954: honor the configured ground plane,
470                // not a hardcoded -1000.0)
471                if state[1] < params.ground_threshold {
472                    break;
473                }
474            }
475
476            // Warn if we hit the iteration limit
477            if iteration_count >= max_iterations {
478                eprintln!(
479                    "WARNING: Trajectory integration hit maximum iteration limit ({} iterations)",
480                    max_iterations
481                );
482                eprintln!("  Final time: {}, Target time: {}", t, t_end);
483                eprintln!(
484                    "  Final position: downrange(x)={}, Target: {}m",
485                    state[0], params.target_distance_m
486                );
487            }
488        }
489    }
490
491    trajectory
492}
493
494/// Python-exposed function for complete trajectory integration
495pub fn solve_trajectory_rust(
496    initial_state: [f64; 6],
497    t_span: (f64, f64),
498    mass_kg: f64,
499    bc: f64,
500    drag_model: DragModel,
501    wind_segments: Vec<WindSegment>,
502    atmos_params: (f64, f64, f64, f64),
503    omega_vector: Option<Vec<f64>>,
504    enable_spin_drift: bool,
505    enable_magnus: bool,
506    enable_coriolis: bool,
507    method: String,
508    tolerance: f64,
509    max_step: f64,
510    target_distance_m: f64,
511) -> Vec<HashMap<String, f64>> {
512    let omega_vec = omega_vector.map(|v| Vector3::new(v[0], v[1], v[2]));
513
514    let params = TrajectoryParams {
515        mass_kg,
516        bc,
517        drag_model,
518        wind_segments,
519        atmos_params,
520        omega_vector: omega_vec,
521        enable_spin_drift,
522        enable_magnus,
523        enable_coriolis,
524        target_distance_m,
525        enable_wind_shear: false, // Default for test function
526        wind_shear_model: "none".to_string(),
527        shooter_altitude_m: 0.0,
528        is_twist_right: true,    // Default for test function
529        custom_drag_table: None, // No CDM for test function
530        bc_segments: None,       // No BC segments for legacy function
531        use_bc_segments: false,
532        ground_threshold: -1000.0, // MBA-954: preserve the historical default
533    };
534
535    let trajectory =
536        integrate_trajectory(initial_state, t_span, params, &method, tolerance, max_step);
537
538    // Convert to Python-friendly format
539    trajectory
540        .into_iter()
541        .map(|(t, state)| {
542            let mut point = HashMap::new();
543            point.insert("t".to_string(), t);
544            point.insert("x".to_string(), state[0]);
545            point.insert("y".to_string(), state[1]);
546            point.insert("z".to_string(), state[2]);
547            point.insert("vx".to_string(), state[3]);
548            point.insert("vy".to_string(), state[4]);
549            point.insert("vz".to_string(), state[5]);
550            point
551        })
552        .collect()
553}
554
555#[cfg(test)]
556mod tests {
557    use super::*;
558
559    fn create_test_params(target_distance_m: f64) -> TrajectoryParams {
560        TrajectoryParams {
561            mass_kg: 0.01134, // 175 grains in kg
562            bc: 0.442,
563            drag_model: DragModel::G7,
564            wind_segments: vec![],
565            atmos_params: (0.0, 59.0, 29.92, 0.0),
566            omega_vector: None,
567            enable_spin_drift: false,
568            enable_magnus: false,
569            enable_coriolis: false,
570            target_distance_m,
571            enable_wind_shear: false,
572            wind_shear_model: "none".to_string(),
573            shooter_altitude_m: 0.0,
574            is_twist_right: true,
575            custom_drag_table: None,
576            bc_segments: None,
577            use_bc_segments: false,
578            ground_threshold: -1000.0,
579        }
580    }
581
582    #[test]
583    fn test_mba954_ground_threshold_honored() {
584        // MBA-954: integrate_trajectory must honor the configured ground plane, not a hardcoded
585        // -1000.0. A descending bullet with a shallow ground_threshold must terminate earlier
586        // (fewer points) than one with the historical deep default.
587        let initial_state = [0.0, 0.0, 0.0, 300.0, -30.0, 0.0]; // descending (vy = -30 m/s)
588
589        let mut shallow = create_test_params(1_000_000.0); // huge target so range never terminates
590        shallow.ground_threshold = -20.0; // stop ~20 m below launch
591        let mut deep = create_test_params(1_000_000.0);
592        deep.ground_threshold = -1000.0; // historical default
593
594        let t_shallow =
595            integrate_trajectory(initial_state, (0.0, 60.0), shallow, "RK4", 1e-6, 0.001);
596        let t_deep = integrate_trajectory(initial_state, (0.0, 60.0), deep, "RK4", 1e-6, 0.001);
597
598        assert!(
599            t_shallow.len() < t_deep.len(),
600            "shallow ground_threshold (-20) should terminate earlier than deep (-1000): \
601             shallow={}, deep={}",
602            t_shallow.len(),
603            t_deep.len()
604        );
605    }
606
607    #[test]
608    fn test_integrate_trajectory_basic() {
609        // Initial state [x,y,z,vx,vy,vz] (McCoy: X=downrange, Z=lateral)
610        // x=0 (downrange start), vx=821.52 (downrange velocity)
611        let initial_state = [0.0, -0.038, 0.0, 821.52, 48.61, 0.0];
612
613        let params = TrajectoryParams {
614            mass_kg: 0.01134, // 175 grains in kg
615            bc: 0.442,
616            drag_model: DragModel::G7,
617            wind_segments: vec![(0.0, 90.0, 914.4)],
618            atmos_params: (0.0, 59.0, 29.92, 0.0),
619            omega_vector: None,
620            enable_spin_drift: false,
621            enable_magnus: false,
622            enable_coriolis: false,
623            target_distance_m: 914.4, // 1000 yards in meters
624            enable_wind_shear: false,
625            wind_shear_model: "none".to_string(),
626            shooter_altitude_m: 0.0,
627            is_twist_right: true,
628            custom_drag_table: None,
629            bc_segments: None,
630            use_bc_segments: false,
631            ground_threshold: -1000.0,
632        };
633
634        println!("Running integrate_trajectory test...");
635        println!("Initial state: {:?}", initial_state);
636        println!("Target distance: {} m", params.target_distance_m);
637
638        let trajectory =
639            integrate_trajectory(initial_state, (0.0, 10.0), params, "RK45", 1e-6, 0.01);
640
641        println!("Trajectory has {} points", trajectory.len());
642
643        // Should have more than just initial point
644        assert!(
645            trajectory.len() > 1,
646            "Trajectory should have more than 1 point, but has {}",
647            trajectory.len()
648        );
649
650        // Check that we actually moved downrange
651        if let Some((_, final_state)) = trajectory.last() {
652            println!("Final state: downrange(x)={}", final_state[0]);
653            assert!(
654                final_state[0] > 0.0,
655                "Final x should be positive (bullet moved downrange)"
656            );
657            assert!(
658                final_state[0] >= 900.0,
659                "Final x should be near target distance"
660            );
661        }
662    }
663
664    #[test]
665    fn test_rk4_vs_rk45_consistency() {
666        // Both methods should give similar results for the same trajectory
667        let initial_state = [0.0, 0.0, 0.0, 800.0, 30.0, 0.0]; // McCoy: vx=downrange
668        let target_distance = 500.0;
669
670        let params_rk4 = create_test_params(target_distance);
671        let params_rk45 = create_test_params(target_distance);
672
673        let trajectory_rk4 =
674            integrate_trajectory(initial_state, (0.0, 5.0), params_rk4, "RK4", 1e-6, 0.001);
675        let trajectory_rk45 =
676            integrate_trajectory(initial_state, (0.0, 5.0), params_rk45, "RK45", 1e-6, 0.01);
677
678        // Both should reach target
679        assert!(!trajectory_rk4.is_empty());
680        assert!(!trajectory_rk45.is_empty());
681
682        let (_, final_rk4) = trajectory_rk4.last().unwrap();
683        let (_, final_rk45) = trajectory_rk45.last().unwrap();
684
685        // Final downrange positions should be within 1% of each other
686        let rk4_z = final_rk4[0];
687        let rk45_z = final_rk45[0];
688        let diff_percent = ((rk4_z - rk45_z) / rk45_z).abs() * 100.0;
689
690        assert!(
691            diff_percent < 1.0,
692            "RK4 and RK45 final positions differ by {}%: RK4={}, RK45={}",
693            diff_percent,
694            rk4_z,
695            rk45_z
696        );
697    }
698
699    #[test]
700    fn test_ground_impact_detection() {
701        // Trajectory with steep downward angle should hit ground
702        let initial_state = [0.0, 100.0, 0.0, 300.0, -50.0, 0.0]; // McCoy: vx=downrange // Steep descent
703
704        let mut params = create_test_params(10000.0); // Far target
705        params.target_distance_m = 10000.0;
706
707        let trajectory =
708            integrate_trajectory(initial_state, (0.0, 20.0), params, "RK45", 1e-6, 0.01);
709
710        // Should stop before reaching target due to ground impact
711        let (_, final_state) = trajectory.last().unwrap();
712
713        // y should be near ground threshold (-1000m)
714        assert!(
715            final_state[1] <= -900.0,
716            "Should hit ground, but y={}",
717            final_state[1]
718        );
719        assert!(
720            final_state[0] < 10000.0,
721            "Should not reach target, but z={}",
722            final_state[0]
723        );
724    }
725
726    #[test]
727    fn test_target_distance_reached() {
728        let initial_state = [0.0, 0.0, 0.0, 800.0, 20.0, 0.0]; // McCoy: vx=downrange
729        let target_distance = 300.0;
730
731        let params = create_test_params(target_distance);
732
733        let trajectory =
734            integrate_trajectory(initial_state, (0.0, 5.0), params, "RK45", 1e-6, 0.01);
735
736        let (_, final_state) = trajectory.last().unwrap();
737
738        // Should stop at or very near target distance
739        assert!(
740            (final_state[0] - target_distance).abs() < 1.0,
741            "Should reach target at {}m, but stopped at {}m",
742            target_distance,
743            final_state[0]
744        );
745    }
746
747    #[test]
748    fn test_wind_affects_trajectory() {
749        // Test that wind segments are properly stored and passed through
750        // The actual wind effect depends on the derivatives computation which
751        // uses the wind vector in the drag calculation
752        let initial_state = [0.0, 0.0, 0.0, 800.0, 30.0, 0.0]; // McCoy: vx=downrange
753        let target_distance = 500.0;
754
755        // No wind
756        let params_no_wind = create_test_params(target_distance);
757
758        // Strong headwind (0 degrees = headwind)
759        let mut params_headwind = create_test_params(target_distance);
760        params_headwind.wind_segments = vec![(72.0, 0.0, 500.0)]; // 72 km/h = 20 m/s headwind
761
762        let trajectory_no_wind =
763            integrate_trajectory(initial_state, (0.0, 5.0), params_no_wind, "RK45", 1e-6, 0.01);
764        let trajectory_headwind =
765            integrate_trajectory(initial_state, (0.0, 5.0), params_headwind, "RK45", 1e-6, 0.01);
766
767        // Both trajectories should complete
768        assert!(!trajectory_no_wind.is_empty(), "No-wind trajectory should complete");
769        assert!(!trajectory_headwind.is_empty(), "Headwind trajectory should complete");
770
771        let (time_no_wind, final_no_wind) = trajectory_no_wind.last().unwrap();
772        let (time_headwind, final_headwind) = trajectory_headwind.last().unwrap();
773
774        // Headwind should slow the bullet, resulting in longer flight time
775        // or different drop at same distance
776        let drop_no_wind = final_no_wind[1];
777        let drop_headwind = final_headwind[1];
778
779        // With headwind, bullet should drop more (more negative y) due to slower velocity
780        // The effect might be small, so we check that both trajectories are valid
781        println!("No wind: time={}, drop={}", time_no_wind, drop_no_wind);
782        println!("Headwind: time={}, drop={}", time_headwind, drop_headwind);
783
784        // Both should reach approximately the target distance
785        assert!(
786            (final_no_wind[0] - target_distance).abs() < 10.0,
787            "No-wind should reach target"
788        );
789        assert!(
790            (final_headwind[0] - target_distance).abs() < 10.0,
791            "Headwind should reach target"
792        );
793    }
794
795    #[test]
796    fn test_solve_trajectory_rust_output_format() {
797        let initial_state = [0.0, 0.0, 0.0, 800.0, 30.0, 0.0]; // McCoy: vx=downrange
798
799        let result = solve_trajectory_rust(
800            initial_state,
801            (0.0, 2.0),
802            0.01134,        // mass_kg
803            0.442,          // bc
804            DragModel::G7,  // drag_model
805            vec![],         // wind_segments
806            (0.0, 59.0, 29.92, 0.0), // atmos_params
807            None,           // omega_vector
808            false,          // enable_spin_drift
809            false,          // enable_magnus
810            false,          // enable_coriolis
811            "RK45".to_string(), // method
812            1e-6,           // tolerance
813            0.01,           // max_step
814            500.0,          // target_distance_m
815        );
816
817        // Should return Vec of HashMaps with expected keys
818        assert!(!result.is_empty());
819
820        let first_point = &result[0];
821        assert!(first_point.contains_key("t"));
822        assert!(first_point.contains_key("x"));
823        assert!(first_point.contains_key("y"));
824        assert!(first_point.contains_key("z"));
825        assert!(first_point.contains_key("vx"));
826        assert!(first_point.contains_key("vy"));
827        assert!(first_point.contains_key("vz"));
828    }
829
830    #[test]
831    fn test_left_vs_right_twist() {
832        let initial_state = [0.0, 0.0, 0.0, 800.0, 30.0, 0.0]; // McCoy: vx=downrange
833        let target_distance = 500.0;
834
835        let mut params_right = create_test_params(target_distance);
836        params_right.is_twist_right = true;
837        params_right.enable_spin_drift = true;
838
839        let mut params_left = create_test_params(target_distance);
840        params_left.is_twist_right = false;
841        params_left.enable_spin_drift = true;
842
843        let trajectory_right =
844            integrate_trajectory(initial_state, (0.0, 5.0), params_right, "RK45", 1e-6, 0.01);
845        let trajectory_left =
846            integrate_trajectory(initial_state, (0.0, 5.0), params_left, "RK45", 1e-6, 0.01);
847
848        // Both should complete
849        assert!(!trajectory_right.is_empty());
850        assert!(!trajectory_left.is_empty());
851
852        // Right and left twist should produce valid trajectories
853        let (_, final_right) = trajectory_right.last().unwrap();
854        let (_, final_left) = trajectory_left.last().unwrap();
855
856        // Both should reach approximately the same downrange distance
857        assert!((final_right[2] - final_left[2]).abs() < 10.0);
858    }
859}