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}
151
152/// Build the loop-invariant BallisticInputs for the derivatives function ONCE per integration,
153/// instead of rebuilding it (a "none".to_string() alloc plus bc_segments / custom_drag_table
154/// clones) on every derivative evaluation (4x per RK4 step, 7x per RK45 step). muzzle_velocity is
155/// set to 0.0 because compute_derivatives never reads it on this path; every other field depends
156/// only on `params`, so the struct is constant for the whole integration.
157fn build_inputs(params: &TrajectoryParams) -> BallisticInputs {
158    BallisticInputs {
159        bc_value: params.bc,
160        bc_type: params.drag_model,
161        bullet_mass: params.mass_kg, // kg
162        muzzle_velocity: 0.0,        // unread by compute_derivatives on this path
163        bullet_diameter: 0.0078232,  // 0.308 in -> meters
164        bullet_length: 0.031496,     // 1.24 in -> meters
165        twist_rate: 10.0,            // inches/turn (twist stays imperial)
166        is_twist_right: params.is_twist_right,
167        enable_advanced_effects: params.enable_spin_drift
168            || params.enable_magnus
169            || params.enable_coriolis,
170        enable_magnus: params.enable_magnus,
171        enable_coriolis: params.enable_coriolis,
172        altitude: params.atmos_params.0,
173        temperature: params.atmos_params.1,
174        pressure: params.atmos_params.2,
175        humidity: params.atmos_params.3,
176        tipoff_yaw: 0.0,
177        target_distance: 1000.0, // default
178        muzzle_angle: 0.0,
179        wind_speed: if !params.wind_segments.is_empty() {
180            params.wind_segments[0].0 * 0.2777778 // km/h -> m/s
181        } else {
182            0.0
183        },
184        wind_angle: if !params.wind_segments.is_empty() {
185            params.wind_segments[0].1.to_radians() // degrees -> radians
186        } else {
187            0.0
188        },
189        latitude: None,
190        shooting_angle: 0.0,
191        azimuth_angle: 0.0,
192        use_powder_sensitivity: false,
193        powder_temp_sensitivity: 0.0,
194        powder_temp: 59.0,
195        tipoff_decay_distance: 0.0,
196        ground_threshold: -1000.0,
197        bc_segments: params.bc_segments.clone(),
198        caliber_inches: 0.308,
199        weight_grains: params.mass_kg / 0.00006479891,
200        use_bc_segments: params.use_bc_segments,
201        bullet_id: None,
202        bc_segments_data: None,
203        use_enhanced_spin_drift: params.enable_spin_drift,
204        use_form_factor: false,
205        manufacturer: None,
206        bullet_model: None,
207        enable_wind_shear: false,
208        wind_shear_model: "none".to_string(),
209        use_cluster_bc: false,
210        bullet_cluster: None,
211        custom_drag_table: params.custom_drag_table.clone(),
212        bc_type_str: None,
213        enable_pitch_damping: false,
214        enable_precession_nutation: false,
215        use_rk4: true,
216        use_adaptive_rk45: false,
217        enable_trajectory_sampling: false,
218        sample_interval: 10.0,
219        sight_height: 0.0,
220        muzzle_height: 0.0,
221        target_height: 0.0,
222    }
223}
224
225/// Convert state to Vector6 and call compute_derivatives
226fn compute_derivatives_vec(
227    state: &Vector6<f64>,
228    t: f64,
229    params: &TrajectoryParams,
230    inputs: &BallisticInputs,
231) -> Vector6<f64> {
232    let pos = Vector3::new(state[0], state[1], state[2]);
233    let vel = Vector3::new(state[3], state[4], state[5]);
234
235    // Calculate wind at current position with shear support
236    let wind_vector = if !params.wind_segments.is_empty() {
237        if params.enable_wind_shear && params.wind_shear_model != "none" {
238            crate::wind_shear::get_wind_at_position(
239                &pos,
240                &params.wind_segments,
241                params.enable_wind_shear,
242                &params.wind_shear_model,
243                params.shooter_altitude_m,
244            )
245        } else {
246            // Simple constant wind (original implementation)
247            let seg = &params.wind_segments[0];
248            let wind_speed_mps = seg.0 * 0.2777778; // km/h to m/s
249            let wind_angle_rad = seg.1.to_radians();
250            // McCoy: x=downrange, y=vertical, z=lateral
251            Vector3::new(
252                -wind_speed_mps * wind_angle_rad.cos(), // x (downrange - head/tail component)
253                0.0,                                    // y (vertical)
254                -wind_speed_mps * wind_angle_rad.sin(), // z (lateral - crosswind component)
255            )
256        }
257    } else {
258        Vector3::zeros()
259    };
260
261    // Call compute_derivatives - returns [f64; 6] directly. `inputs` is built once per
262    // integration by build_inputs() and threaded in, instead of rebuilt every call.
263    let deriv_result = compute_derivatives(
264        pos,
265        vel,
266        inputs,
267        wind_vector,
268        params.atmos_params,
269        params.bc,
270        params.omega_vector,
271        t,
272    );
273
274    Vector6::new(
275        deriv_result[0],
276        deriv_result[1],
277        deriv_result[2],
278        deriv_result[3],
279        deriv_result[4],
280        deriv_result[5],
281    )
282}
283
284/// Main trajectory integration function
285pub fn integrate_trajectory(
286    initial_state: [f64; 6],
287    t_span: (f64, f64),
288    params: TrajectoryParams,
289    method: &str,
290    tolerance: f64,
291    max_step: f64,
292) -> Vec<(f64, Vector6<f64>)> {
293    let mut state = Vector6::new(
294        initial_state[0],
295        initial_state[1],
296        initial_state[2],
297        initial_state[3],
298        initial_state[4],
299        initial_state[5],
300    );
301
302    let mut t = t_span.0;
303    let t_end = t_span.1;
304    let mut dt = (t_end - t) / 1000.0; // Initial step size
305
306    let mut trajectory = Vec::with_capacity(10000);
307    trajectory.push((t, state));
308
309    // Build the (loop-invariant) derivative inputs once for the whole integration, instead of
310    // rebuilding the struct on every derivative evaluation.
311    let inputs = build_inputs(&params);
312
313    match method {
314        "RK4" => {
315            // Fixed step RK4 with target detection
316            dt = dt.min(max_step).min(0.001); // Use smaller steps for accuracy
317
318            while t < t_end {
319                if t + dt > t_end {
320                    dt = t_end - t;
321                }
322
323                let new_state = rk4_step(&state, t, dt, &params, &inputs);
324
325                // Check if we're about to pass the target (X is downrange, McCoy)
326                if state[0] < params.target_distance_m && new_state[0] >= params.target_distance_m {
327                    // Interpolate to find exact target crossing
328                    let alpha = (params.target_distance_m - state[0]) / (new_state[0] - state[0]);
329                    let dt_to_target = dt * alpha;
330
331                    // Take a smaller step to reach target exactly
332                    let final_state = rk4_step(&state, t, dt_to_target, &params, &inputs);
333
334                    // Ensure we don't overshoot
335                    let mut corrected_state = final_state;
336                    if corrected_state[0] > params.target_distance_m {
337                        corrected_state[0] = params.target_distance_m;
338                    }
339
340                    trajectory.push((t + dt_to_target, corrected_state));
341                    break; // Stop at target
342                }
343
344                state = new_state;
345                t += dt;
346                trajectory.push((t, state));
347
348                // Check if we've reached or passed the target
349                if state[0] >= params.target_distance_m {
350                    // X is downrange (McCoy)
351                    // Add final point exactly at target
352                    let mut final_state = state;
353                    final_state[0] = params.target_distance_m; // X is downrange
354                    trajectory.push((t, final_state));
355                    break;
356                }
357
358                // Check if bullet hit ground
359                if state[1] < -1000.0 {
360                    break;
361                }
362            }
363        }
364        "RK45" | _ => {
365            // Adaptive RK45 with better sampling
366            let mut last_save_x = 0.0; // X is downrange (McCoy)
367            let save_interval_m = params.target_distance_m / 50.0; // Save ~50 points minimum
368
369            // OPTIMIZATION: Adjust max step size when wind shear is enabled
370            // This improves numerical stability at long ranges
371            let effective_max_step =
372                if params.enable_wind_shear && params.wind_shear_model != "none" {
373                    // Use smaller steps for wind shear, but not TOO small
374                    if params.target_distance_m > 800.0 {
375                        0.01 // Smaller steps for long range with shear (10ms)
376                    } else {
377                        0.02 // Normal steps for medium range with shear (20ms)
378                    }
379                } else {
380                    max_step // Use provided max_step when no wind shear
381                };
382
383            // Set initial step size - ensure it's reasonable
384            dt = dt.min(effective_max_step).max(0.0001); // At least 0.1ms to avoid infinite loops
385
386            // Safety check: maximum iterations to prevent infinite loops
387            let max_iterations = 100000; // Should be more than enough for any realistic trajectory
388            let mut iteration_count = 0;
389
390            while t < t_end && iteration_count < max_iterations {
391                iteration_count += 1;
392
393                // Limit time step for better resolution
394                if t + dt > t_end {
395                    dt = t_end - t;
396                }
397
398                let (new_state, dt_new, _error) =
399                    rk45_step(&state, t, dt, &params, &inputs, tolerance);
400
401                // Check if we're about to pass the target (X is downrange, McCoy)
402                if state[0] < params.target_distance_m && new_state[0] >= params.target_distance_m {
403                    // Interpolate to find exact target crossing
404                    let alpha = (params.target_distance_m - state[0]) / (new_state[0] - state[0]);
405                    let dt_to_target = dt * alpha;
406
407                    // Take a smaller step to reach target exactly
408                    let (final_state, _, _) =
409                        rk45_step(&state, t, dt_to_target, &params, &inputs, tolerance);
410
411                    // Make sure we don't overshoot
412                    let mut corrected_state = final_state;
413                    if corrected_state[0] > params.target_distance_m {
414                        corrected_state[0] = params.target_distance_m;
415                    }
416
417                    trajectory.push((t + dt_to_target, corrected_state));
418                    break; // Stop at target - no more points after this
419                }
420
421                // Update state
422                state = new_state;
423                t += dt;
424
425                // Save trajectory point if we've moved enough distance
426                if state[0] - last_save_x >= save_interval_m || state[0] >= params.target_distance_m
427                {
428                    // X is downrange
429                    trajectory.push((t, state));
430                    last_save_x = state[0];
431                }
432
433                // Limit dt for next step - ensure we get enough resolution
434                dt = dt_new.min(effective_max_step).max(0.0001); // Use effective max step, min 0.1ms
435
436                // Stop if we've reached the target
437                if state[0] >= params.target_distance_m {
438                    // X is downrange (McCoy)
439                    // Add final point at target distance
440                    let mut final_state = state;
441                    final_state[0] = params.target_distance_m; // X is downrange
442                    trajectory.push((t, final_state));
443                    break;
444                }
445
446                // Check if bullet hit ground
447                if state[1] < -1000.0 {
448                    break;
449                }
450            }
451
452            // Warn if we hit the iteration limit
453            if iteration_count >= max_iterations {
454                eprintln!(
455                    "WARNING: Trajectory integration hit maximum iteration limit ({} iterations)",
456                    max_iterations
457                );
458                eprintln!("  Final time: {}, Target time: {}", t, t_end);
459                eprintln!(
460                    "  Final position: downrange(x)={}, Target: {}m",
461                    state[0], params.target_distance_m
462                );
463            }
464        }
465    }
466
467    trajectory
468}
469
470/// Python-exposed function for complete trajectory integration
471pub fn solve_trajectory_rust(
472    initial_state: [f64; 6],
473    t_span: (f64, f64),
474    mass_kg: f64,
475    bc: f64,
476    drag_model: DragModel,
477    wind_segments: Vec<WindSegment>,
478    atmos_params: (f64, f64, f64, f64),
479    omega_vector: Option<Vec<f64>>,
480    enable_spin_drift: bool,
481    enable_magnus: bool,
482    enable_coriolis: bool,
483    method: String,
484    tolerance: f64,
485    max_step: f64,
486    target_distance_m: f64,
487) -> Vec<HashMap<String, f64>> {
488    let omega_vec = omega_vector.map(|v| Vector3::new(v[0], v[1], v[2]));
489
490    let params = TrajectoryParams {
491        mass_kg,
492        bc,
493        drag_model,
494        wind_segments,
495        atmos_params,
496        omega_vector: omega_vec,
497        enable_spin_drift,
498        enable_magnus,
499        enable_coriolis,
500        target_distance_m,
501        enable_wind_shear: false, // Default for test function
502        wind_shear_model: "none".to_string(),
503        shooter_altitude_m: 0.0,
504        is_twist_right: true,    // Default for test function
505        custom_drag_table: None, // No CDM for test function
506        bc_segments: None,       // No BC segments for legacy function
507        use_bc_segments: false,
508    };
509
510    let trajectory =
511        integrate_trajectory(initial_state, t_span, params, &method, tolerance, max_step);
512
513    // Convert to Python-friendly format
514    trajectory
515        .into_iter()
516        .map(|(t, state)| {
517            let mut point = HashMap::new();
518            point.insert("t".to_string(), t);
519            point.insert("x".to_string(), state[0]);
520            point.insert("y".to_string(), state[1]);
521            point.insert("z".to_string(), state[2]);
522            point.insert("vx".to_string(), state[3]);
523            point.insert("vy".to_string(), state[4]);
524            point.insert("vz".to_string(), state[5]);
525            point
526        })
527        .collect()
528}
529
530#[cfg(test)]
531mod tests {
532    use super::*;
533
534    fn create_test_params(target_distance_m: f64) -> TrajectoryParams {
535        TrajectoryParams {
536            mass_kg: 0.01134, // 175 grains in kg
537            bc: 0.442,
538            drag_model: DragModel::G7,
539            wind_segments: vec![],
540            atmos_params: (0.0, 59.0, 29.92, 0.0),
541            omega_vector: None,
542            enable_spin_drift: false,
543            enable_magnus: false,
544            enable_coriolis: false,
545            target_distance_m,
546            enable_wind_shear: false,
547            wind_shear_model: "none".to_string(),
548            shooter_altitude_m: 0.0,
549            is_twist_right: true,
550            custom_drag_table: None,
551            bc_segments: None,
552            use_bc_segments: false,
553        }
554    }
555
556    #[test]
557    fn test_integrate_trajectory_basic() {
558        // Initial state [x,y,z,vx,vy,vz] (McCoy: X=downrange, Z=lateral)
559        // x=0 (downrange start), vx=821.52 (downrange velocity)
560        let initial_state = [0.0, -0.038, 0.0, 821.52, 48.61, 0.0];
561
562        let params = TrajectoryParams {
563            mass_kg: 0.01134, // 175 grains in kg
564            bc: 0.442,
565            drag_model: DragModel::G7,
566            wind_segments: vec![(0.0, 90.0, 914.4)],
567            atmos_params: (0.0, 59.0, 29.92, 0.0),
568            omega_vector: None,
569            enable_spin_drift: false,
570            enable_magnus: false,
571            enable_coriolis: false,
572            target_distance_m: 914.4, // 1000 yards in meters
573            enable_wind_shear: false,
574            wind_shear_model: "none".to_string(),
575            shooter_altitude_m: 0.0,
576            is_twist_right: true,
577            custom_drag_table: None,
578            bc_segments: None,
579            use_bc_segments: false,
580        };
581
582        println!("Running integrate_trajectory test...");
583        println!("Initial state: {:?}", initial_state);
584        println!("Target distance: {} m", params.target_distance_m);
585
586        let trajectory =
587            integrate_trajectory(initial_state, (0.0, 10.0), params, "RK45", 1e-6, 0.01);
588
589        println!("Trajectory has {} points", trajectory.len());
590
591        // Should have more than just initial point
592        assert!(
593            trajectory.len() > 1,
594            "Trajectory should have more than 1 point, but has {}",
595            trajectory.len()
596        );
597
598        // Check that we actually moved downrange
599        if let Some((_, final_state)) = trajectory.last() {
600            println!("Final state: downrange(x)={}", final_state[0]);
601            assert!(
602                final_state[0] > 0.0,
603                "Final x should be positive (bullet moved downrange)"
604            );
605            assert!(
606                final_state[0] >= 900.0,
607                "Final x should be near target distance"
608            );
609        }
610    }
611
612    #[test]
613    fn test_rk4_vs_rk45_consistency() {
614        // Both methods should give similar results for the same trajectory
615        let initial_state = [0.0, 0.0, 0.0, 800.0, 30.0, 0.0]; // McCoy: vx=downrange
616        let target_distance = 500.0;
617
618        let params_rk4 = create_test_params(target_distance);
619        let params_rk45 = create_test_params(target_distance);
620
621        let trajectory_rk4 =
622            integrate_trajectory(initial_state, (0.0, 5.0), params_rk4, "RK4", 1e-6, 0.001);
623        let trajectory_rk45 =
624            integrate_trajectory(initial_state, (0.0, 5.0), params_rk45, "RK45", 1e-6, 0.01);
625
626        // Both should reach target
627        assert!(!trajectory_rk4.is_empty());
628        assert!(!trajectory_rk45.is_empty());
629
630        let (_, final_rk4) = trajectory_rk4.last().unwrap();
631        let (_, final_rk45) = trajectory_rk45.last().unwrap();
632
633        // Final downrange positions should be within 1% of each other
634        let rk4_z = final_rk4[0];
635        let rk45_z = final_rk45[0];
636        let diff_percent = ((rk4_z - rk45_z) / rk45_z).abs() * 100.0;
637
638        assert!(
639            diff_percent < 1.0,
640            "RK4 and RK45 final positions differ by {}%: RK4={}, RK45={}",
641            diff_percent,
642            rk4_z,
643            rk45_z
644        );
645    }
646
647    #[test]
648    fn test_ground_impact_detection() {
649        // Trajectory with steep downward angle should hit ground
650        let initial_state = [0.0, 100.0, 0.0, 300.0, -50.0, 0.0]; // McCoy: vx=downrange // Steep descent
651
652        let mut params = create_test_params(10000.0); // Far target
653        params.target_distance_m = 10000.0;
654
655        let trajectory =
656            integrate_trajectory(initial_state, (0.0, 20.0), params, "RK45", 1e-6, 0.01);
657
658        // Should stop before reaching target due to ground impact
659        let (_, final_state) = trajectory.last().unwrap();
660
661        // y should be near ground threshold (-1000m)
662        assert!(
663            final_state[1] <= -900.0,
664            "Should hit ground, but y={}",
665            final_state[1]
666        );
667        assert!(
668            final_state[0] < 10000.0,
669            "Should not reach target, but z={}",
670            final_state[0]
671        );
672    }
673
674    #[test]
675    fn test_target_distance_reached() {
676        let initial_state = [0.0, 0.0, 0.0, 800.0, 20.0, 0.0]; // McCoy: vx=downrange
677        let target_distance = 300.0;
678
679        let params = create_test_params(target_distance);
680
681        let trajectory =
682            integrate_trajectory(initial_state, (0.0, 5.0), params, "RK45", 1e-6, 0.01);
683
684        let (_, final_state) = trajectory.last().unwrap();
685
686        // Should stop at or very near target distance
687        assert!(
688            (final_state[0] - target_distance).abs() < 1.0,
689            "Should reach target at {}m, but stopped at {}m",
690            target_distance,
691            final_state[0]
692        );
693    }
694
695    #[test]
696    fn test_wind_affects_trajectory() {
697        // Test that wind segments are properly stored and passed through
698        // The actual wind effect depends on the derivatives computation which
699        // uses the wind vector in the drag calculation
700        let initial_state = [0.0, 0.0, 0.0, 800.0, 30.0, 0.0]; // McCoy: vx=downrange
701        let target_distance = 500.0;
702
703        // No wind
704        let params_no_wind = create_test_params(target_distance);
705
706        // Strong headwind (0 degrees = headwind)
707        let mut params_headwind = create_test_params(target_distance);
708        params_headwind.wind_segments = vec![(72.0, 0.0, 500.0)]; // 72 km/h = 20 m/s headwind
709
710        let trajectory_no_wind =
711            integrate_trajectory(initial_state, (0.0, 5.0), params_no_wind, "RK45", 1e-6, 0.01);
712        let trajectory_headwind =
713            integrate_trajectory(initial_state, (0.0, 5.0), params_headwind, "RK45", 1e-6, 0.01);
714
715        // Both trajectories should complete
716        assert!(!trajectory_no_wind.is_empty(), "No-wind trajectory should complete");
717        assert!(!trajectory_headwind.is_empty(), "Headwind trajectory should complete");
718
719        let (time_no_wind, final_no_wind) = trajectory_no_wind.last().unwrap();
720        let (time_headwind, final_headwind) = trajectory_headwind.last().unwrap();
721
722        // Headwind should slow the bullet, resulting in longer flight time
723        // or different drop at same distance
724        let drop_no_wind = final_no_wind[1];
725        let drop_headwind = final_headwind[1];
726
727        // With headwind, bullet should drop more (more negative y) due to slower velocity
728        // The effect might be small, so we check that both trajectories are valid
729        println!("No wind: time={}, drop={}", time_no_wind, drop_no_wind);
730        println!("Headwind: time={}, drop={}", time_headwind, drop_headwind);
731
732        // Both should reach approximately the target distance
733        assert!(
734            (final_no_wind[0] - target_distance).abs() < 10.0,
735            "No-wind should reach target"
736        );
737        assert!(
738            (final_headwind[0] - target_distance).abs() < 10.0,
739            "Headwind should reach target"
740        );
741    }
742
743    #[test]
744    fn test_solve_trajectory_rust_output_format() {
745        let initial_state = [0.0, 0.0, 0.0, 800.0, 30.0, 0.0]; // McCoy: vx=downrange
746
747        let result = solve_trajectory_rust(
748            initial_state,
749            (0.0, 2.0),
750            0.01134,        // mass_kg
751            0.442,          // bc
752            DragModel::G7,  // drag_model
753            vec![],         // wind_segments
754            (0.0, 59.0, 29.92, 0.0), // atmos_params
755            None,           // omega_vector
756            false,          // enable_spin_drift
757            false,          // enable_magnus
758            false,          // enable_coriolis
759            "RK45".to_string(), // method
760            1e-6,           // tolerance
761            0.01,           // max_step
762            500.0,          // target_distance_m
763        );
764
765        // Should return Vec of HashMaps with expected keys
766        assert!(!result.is_empty());
767
768        let first_point = &result[0];
769        assert!(first_point.contains_key("t"));
770        assert!(first_point.contains_key("x"));
771        assert!(first_point.contains_key("y"));
772        assert!(first_point.contains_key("z"));
773        assert!(first_point.contains_key("vx"));
774        assert!(first_point.contains_key("vy"));
775        assert!(first_point.contains_key("vz"));
776    }
777
778    #[test]
779    fn test_left_vs_right_twist() {
780        let initial_state = [0.0, 0.0, 0.0, 800.0, 30.0, 0.0]; // McCoy: vx=downrange
781        let target_distance = 500.0;
782
783        let mut params_right = create_test_params(target_distance);
784        params_right.is_twist_right = true;
785        params_right.enable_spin_drift = true;
786
787        let mut params_left = create_test_params(target_distance);
788        params_left.is_twist_right = false;
789        params_left.enable_spin_drift = true;
790
791        let trajectory_right =
792            integrate_trajectory(initial_state, (0.0, 5.0), params_right, "RK45", 1e-6, 0.01);
793        let trajectory_left =
794            integrate_trajectory(initial_state, (0.0, 5.0), params_left, "RK45", 1e-6, 0.01);
795
796        // Both should complete
797        assert!(!trajectory_right.is_empty());
798        assert!(!trajectory_left.is_empty());
799
800        // Right and left twist should produce valid trajectories
801        let (_, final_right) = trajectory_right.last().unwrap();
802        let (_, final_left) = trajectory_left.last().unwrap();
803
804        // Both should reach approximately the same downrange distance
805        assert!((final_right[2] - final_left[2]).abs() < 10.0);
806    }
807}