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