Skip to main content

pavan/
vlm.rs

1//! Vortex Lattice Method (VLM) for 3D finite wings.
2//!
3//! Models a wing as a lattice of horseshoe vortices, solving for circulation
4//! distribution via Biot-Savart influence coefficients. Produces span loading,
5//! lift coefficient, induced drag, and Oswald efficiency factor.
6
7use std::f64::consts::PI;
8
9use hisab::DVec3;
10use serde::{Deserialize, Serialize};
11
12use crate::error::{PavanError, Result};
13
14/// Vortex core cutoff to prevent singularity in Biot-Savart.
15const VORTEX_CORE_EPS: f64 = 1e-8;
16
17/// Far-field truncation length for semi-infinite trailing legs (chord-lengths).
18const TRAILING_LEG_FACTOR: f64 = 50.0;
19
20/// Minimum number of panels for a meaningful solve.
21const MIN_VLM_PANELS: usize = 2;
22
23/// Wing planform definition for VLM analysis.
24#[derive(Debug, Clone, Serialize, Deserialize)]
25#[non_exhaustive]
26pub struct WingGeometry {
27    /// Full wingspan (m), tip to tip.
28    pub span: f64,
29    /// Root chord length (m).
30    pub root_chord: f64,
31    /// Tip chord length (m). Equal to root_chord for rectangular wings.
32    pub tip_chord: f64,
33    /// Leading edge sweep angle (rad). 0 for unswept.
34    pub sweep_le_rad: f64,
35    /// Dihedral angle (rad). 0 for flat wing.
36    pub dihedral_rad: f64,
37    /// Geometric twist at tip relative to root (rad). Negative = washout.
38    pub twist_tip_rad: f64,
39    /// Spanwise panels per semi-span.
40    pub n_span: usize,
41    /// Chordwise panels.
42    pub n_chord: usize,
43}
44
45impl WingGeometry {
46    /// Rectangular wing (no taper, no sweep, no twist).
47    #[must_use]
48    #[inline]
49    pub fn rectangular(span: f64, chord: f64, n_span: usize, n_chord: usize) -> Self {
50        Self {
51            span,
52            root_chord: chord,
53            tip_chord: chord,
54            sweep_le_rad: 0.0,
55            dihedral_rad: 0.0,
56            twist_tip_rad: 0.0,
57            n_span,
58            n_chord,
59        }
60    }
61
62    /// Tapered wing (no sweep, no twist).
63    #[must_use]
64    #[inline]
65    pub fn tapered(
66        span: f64,
67        root_chord: f64,
68        tip_chord: f64,
69        n_span: usize,
70        n_chord: usize,
71    ) -> Self {
72        Self {
73            span,
74            root_chord,
75            tip_chord,
76            sweep_le_rad: 0.0,
77            dihedral_rad: 0.0,
78            twist_tip_rad: 0.0,
79            n_span,
80            n_chord,
81        }
82    }
83
84    /// Wing reference area: S = span × (root + tip) / 2.
85    #[must_use]
86    #[inline]
87    pub fn reference_area(&self) -> f64 {
88        self.span * (self.root_chord + self.tip_chord) * 0.5
89    }
90
91    /// Aspect ratio: AR = span² / S.
92    #[must_use]
93    #[inline]
94    pub fn aspect_ratio(&self) -> f64 {
95        let s = self.reference_area();
96        if s <= 0.0 {
97            return 0.0;
98        }
99        self.span * self.span / s
100    }
101
102    /// Total number of panels (both semi-spans × chordwise).
103    #[must_use]
104    #[inline]
105    pub fn total_panels(&self) -> usize {
106        2 * self.n_span * self.n_chord
107    }
108}
109
110/// A single VLM panel with horseshoe vortex.
111#[derive(Debug, Clone)]
112pub struct VlmPanel {
113    /// Panel corners: [LE-left, LE-right, TE-right, TE-left].
114    pub corners: [DVec3; 4],
115    /// Bound vortex endpoints (at 1/4 chord): [left, right].
116    pub bound: [DVec3; 2],
117    /// Control point at 3/4 chord, mid-span.
118    pub control_point: DVec3,
119    /// Outward unit normal.
120    pub normal: DVec3,
121    /// Panel area (m²).
122    pub area: f64,
123    /// Local chord length (m).
124    pub chord: f64,
125    /// Spanwise width (m).
126    pub dy: f64,
127}
128
129/// Result of a VLM solve.
130#[derive(Debug, Clone, Serialize, Deserialize)]
131#[non_exhaustive]
132pub struct VlmSolution {
133    /// Total lift coefficient.
134    pub cl: f64,
135    /// Induced drag coefficient.
136    pub cdi: f64,
137    /// Pitching moment coefficient about root leading edge.
138    pub cm: f64,
139    /// Sectional lift coefficient per spanwise strip.
140    pub span_cl: Vec<f64>,
141    /// Circulation strength per panel (m²/s).
142    pub circulation: Vec<f64>,
143    /// Oswald span efficiency factor (1.0 = elliptic loading).
144    pub oswald_efficiency: f64,
145}
146
147/// Generate VLM panels from wing planform geometry.
148///
149/// Creates a lattice of panels across both semi-spans, with bound vortices
150/// at 1/4 chord and control points at 3/4 chord (Pistolesi's theorem).
151#[must_use]
152pub fn generate_panels(wing: &WingGeometry) -> Vec<VlmPanel> {
153    let n_span_total = 2 * wing.n_span;
154    let n_chord = wing.n_chord;
155    let n_total = n_span_total * n_chord;
156    let mut panels = Vec::with_capacity(n_total);
157
158    let half_span = wing.span * 0.5;
159
160    for i_span in 0..n_span_total {
161        // Spanwise position: -half_span to +half_span
162        let eta_left = (i_span as f64 / n_span_total as f64) * 2.0 - 1.0;
163        let eta_right = ((i_span + 1) as f64 / n_span_total as f64) * 2.0 - 1.0;
164        let y_left = eta_left * half_span;
165        let y_right = eta_right * half_span;
166
167        // Interpolate chord and twist at each spanwise station
168        let frac_left = eta_left.abs();
169        let frac_right = eta_right.abs();
170        let chord_left = wing.root_chord + (wing.tip_chord - wing.root_chord) * frac_left;
171        let chord_right = wing.root_chord + (wing.tip_chord - wing.root_chord) * frac_right;
172
173        // Leading edge x-offset from sweep
174        let x_le_left = frac_left * half_span * wing.sweep_le_rad.tan();
175        let x_le_right = frac_right * half_span * wing.sweep_le_rad.tan();
176
177        // Dihedral z-offset
178        let z_left = y_left.abs() * wing.dihedral_rad.sin();
179        let z_right = y_right.abs() * wing.dihedral_rad.sin();
180
181        // Twist (linear from root to tip)
182        let twist_left = wing.twist_tip_rad * frac_left;
183        let twist_right = wing.twist_tip_rad * frac_right;
184
185        for i_chord in 0..n_chord {
186            let xi_le = i_chord as f64 / n_chord as f64;
187            let xi_te = (i_chord + 1) as f64 / n_chord as f64;
188
189            // Panel corners (x downstream, y spanwise, z up)
190            let le_left = DVec3::new(
191                x_le_left + xi_le * chord_left,
192                y_left,
193                z_left - xi_le * chord_left * twist_left.sin(),
194            );
195            let le_right = DVec3::new(
196                x_le_right + xi_le * chord_right,
197                y_right,
198                z_right - xi_le * chord_right * twist_right.sin(),
199            );
200            let te_right = DVec3::new(
201                x_le_right + xi_te * chord_right,
202                y_right,
203                z_right - xi_te * chord_right * twist_right.sin(),
204            );
205            let te_left = DVec3::new(
206                x_le_left + xi_te * chord_left,
207                y_left,
208                z_left - xi_te * chord_left * twist_left.sin(),
209            );
210
211            // Bound vortex at 1/4 chord of this panel
212            let quarter = 0.25;
213            let bound_frac = xi_le + quarter * (xi_te - xi_le);
214            let bound_left = DVec3::new(
215                x_le_left + bound_frac * chord_left,
216                y_left,
217                z_left - bound_frac * chord_left * twist_left.sin(),
218            );
219            let bound_right = DVec3::new(
220                x_le_right + bound_frac * chord_right,
221                y_right,
222                z_right - bound_frac * chord_right * twist_right.sin(),
223            );
224
225            // Control point at 3/4 chord, mid-span
226            let three_quarter = 0.75;
227            let cp_frac = xi_le + three_quarter * (xi_te - xi_le);
228            let y_mid = 0.5 * (y_left + y_right);
229            let frac_mid = y_mid.abs() / half_span;
230            let chord_mid = wing.root_chord + (wing.tip_chord - wing.root_chord) * frac_mid;
231            let x_le_mid = frac_mid * half_span * wing.sweep_le_rad.tan();
232            let z_mid = y_mid.abs() * wing.dihedral_rad.sin();
233            let twist_mid = wing.twist_tip_rad * frac_mid;
234            let cp = DVec3::new(
235                x_le_mid + cp_frac * chord_mid,
236                y_mid,
237                z_mid - cp_frac * chord_mid * twist_mid.sin(),
238            );
239
240            // Normal from cross product of diagonals (order gives +Z for flat wing)
241            let diag1 = te_right - le_left;
242            let diag2 = te_left - le_right;
243            let normal = diag2.cross(diag1);
244            let normal_len = normal.length();
245            let normal = if normal_len > f64::EPSILON {
246                normal / normal_len
247            } else {
248                tracing::warn!(
249                    "degenerate VLM panel at span index {i_span}, chord index {i_chord}: zero normal"
250                );
251                DVec3::Z
252            };
253
254            // Area = 0.5 * |diag1 × diag2|
255            let area = 0.5 * normal_len;
256            let local_chord = 0.5 * (chord_left + chord_right) * (xi_te - xi_le);
257            let dy = (y_right - y_left).abs();
258
259            panels.push(VlmPanel {
260                corners: [le_left, le_right, te_right, te_left],
261                bound: [bound_left, bound_right],
262                control_point: cp,
263                normal,
264                area,
265                chord: local_chord,
266                dy,
267            });
268        }
269    }
270
271    panels
272}
273
274/// Velocity induced by a finite vortex segment a→b at point p (unit circulation).
275///
276/// Uses the Biot-Savart law with vortex core cutoff for numerical stability.
277#[must_use]
278#[inline]
279pub fn biot_savart_segment(p: DVec3, a: DVec3, b: DVec3) -> DVec3 {
280    let r1 = p - a;
281    let r2 = p - b;
282    let r0 = b - a;
283
284    let cross = r1.cross(r2);
285    let cross_sq = cross.length_squared();
286
287    // Core cutoff: avoid singularity when p is on the vortex line
288    if cross_sq < VORTEX_CORE_EPS {
289        return DVec3::ZERO;
290    }
291
292    let r1_len = r1.length();
293    let r2_len = r2.length();
294
295    if r1_len < VORTEX_CORE_EPS || r2_len < VORTEX_CORE_EPS {
296        return DVec3::ZERO;
297    }
298
299    // Biot-Savart: V = Γ/(4π) · (r1×r2)/|r1×r2|² · r0·(r1/|r1| - r2/|r2|)
300    let dot_factor = r0.dot(r1 / r1_len - r2 / r2_len);
301
302    (1.0 / (4.0 * PI)) * cross * (dot_factor / cross_sq)
303}
304
305/// Induced velocity from a horseshoe vortex (unit circulation) at a control point.
306///
307/// Horseshoe = bound segment + two trailing semi-infinite legs extending downstream (+x).
308#[inline]
309fn horseshoe_influence(panel: &VlmPanel, control: DVec3, far_field: f64) -> DVec3 {
310    let a = panel.bound[0]; // left end of bound vortex
311    let b = panel.bound[1]; // right end of bound vortex
312
313    // Trailing leg direction: downstream (+x)
314    let trail = DVec3::new(far_field, 0.0, 0.0);
315
316    // Left trailing leg: from far upstream-left to bound-left (semi-infinite, approximated)
317    let left_far = a + trail;
318    let v_left = biot_savart_segment(control, left_far, a);
319
320    // Bound segment: from left to right
321    let v_bound = biot_savart_segment(control, a, b);
322
323    // Right trailing leg: from bound-right to far downstream-right
324    let right_far = b + trail;
325    let v_right = biot_savart_segment(control, b, right_far);
326
327    v_left + v_bound + v_right
328}
329
330/// Post-process VLM solution: compute CL, CDi, Cm, span loading, Oswald efficiency.
331#[allow(clippy::needless_range_loop)]
332fn vlm_post_process(
333    panels: &[VlmPanel],
334    gamma: Vec<f64>,
335    wing: &WingGeometry,
336    v_inf: f64,
337) -> VlmSolution {
338    let s_ref = wing.reference_area();
339    let q_inf = 0.5 * v_inf * v_inf;
340    let ar = wing.aspect_ratio();
341    let half_span = wing.span * 0.5;
342    let n_span_total = 2 * wing.n_span;
343    let n_chord = wing.n_chord;
344
345    // Aggregate circulation per spanwise strip
346    let mut strip_gamma = vec![0.0; n_span_total];
347    for i_span in 0..n_span_total {
348        for i_chord in 0..n_chord {
349            strip_gamma[i_span] += gamma[i_span * n_chord + i_chord];
350        }
351    }
352
353    // Lift and moment via Kutta-Joukowski
354    let mut cl_total = 0.0;
355    let mut cm_total = 0.0;
356    let mut span_cl = Vec::with_capacity(n_span_total);
357
358    for i_span in 0..n_span_total {
359        let idx0 = i_span * n_chord;
360        let dy = panels[idx0].dy;
361        let local_chord: f64 = (0..n_chord).map(|ic| panels[idx0 + ic].chord).sum();
362
363        let cl_section = if local_chord > 0.0 && q_inf > 0.0 {
364            strip_gamma[i_span] * v_inf / (q_inf * local_chord)
365        } else {
366            0.0
367        };
368        span_cl.push(cl_section);
369
370        let dcl = strip_gamma[i_span] * dy * v_inf / (q_inf * s_ref);
371        cl_total += dcl;
372
373        let x_cp = panels[idx0].control_point.x;
374        let mac = 0.5 * (wing.root_chord + wing.tip_chord);
375        cm_total -= x_cp * dcl / mac;
376    }
377
378    // Induced drag via Trefftz plane trailing vortex sheet
379    let mut trail_gamma = Vec::with_capacity(n_span_total + 1);
380    trail_gamma.push(strip_gamma[0]);
381    for j in 0..n_span_total - 1 {
382        trail_gamma.push(strip_gamma[j + 1] - strip_gamma[j]);
383    }
384    trail_gamma.push(-strip_gamma[n_span_total - 1]);
385
386    let mut trail_y = Vec::with_capacity(n_span_total + 1);
387    for j in 0..=n_span_total {
388        let eta = (j as f64 / n_span_total as f64) * 2.0 - 1.0;
389        trail_y.push(eta * half_span);
390    }
391
392    let mut cdi_total = 0.0;
393    for i_span in 0..n_span_total {
394        let idx0 = i_span * n_chord;
395        let y_i = panels[idx0].control_point.y;
396        let dy = panels[idx0].dy;
397
398        let mut w_i = 0.0;
399        for j in 0..trail_gamma.len() {
400            let dist = y_i - trail_y[j];
401            if dist.abs() > 1e-10 {
402                w_i -= trail_gamma[j] / (4.0 * PI * dist);
403            }
404        }
405        cdi_total -= strip_gamma[i_span] * w_i * dy / (q_inf * s_ref);
406    }
407
408    let oswald = if cdi_total.abs() > f64::EPSILON && ar > 0.0 {
409        cl_total * cl_total / (PI * ar * cdi_total)
410    } else {
411        1.0
412    };
413
414    VlmSolution {
415        cl: cl_total,
416        cdi: cdi_total,
417        cm: cm_total,
418        span_cl,
419        circulation: gamma,
420        oswald_efficiency: oswald,
421    }
422}
423
424/// Solve the VLM for a single angle of attack.
425///
426/// Panels should be generated from [`generate_panels`]. Freestream velocity
427/// `v_inf` in m/s; angle of attack in radians.
428#[allow(clippy::needless_range_loop)]
429pub fn solve(
430    panels: &[VlmPanel],
431    wing: &WingGeometry,
432    alpha_rad: f64,
433    v_inf: f64,
434) -> Result<VlmSolution> {
435    let n = panels.len();
436    if n < MIN_VLM_PANELS {
437        return Err(PavanError::InvalidGeometry(format!(
438            "need at least {MIN_VLM_PANELS} VLM panels, got {n}"
439        )));
440    }
441    if v_inf <= 0.0 {
442        return Err(PavanError::InvalidVelocity(
443            "freestream velocity must be positive".into(),
444        ));
445    }
446
447    let far_field = wing.root_chord * TRAILING_LEG_FACTOR;
448
449    // Freestream velocity vector
450    let v_free = DVec3::new(v_inf * alpha_rad.cos(), 0.0, v_inf * alpha_rad.sin());
451
452    // Build influence matrix: A[i][j] = V_induced_by_horseshoe_j at control_i · normal_i
453    let mut matrix: Vec<Vec<f64>> = vec![vec![0.0; n]; n];
454    let mut rhs: Vec<f64> = vec![0.0; n];
455
456    for i in 0..n {
457        for j in 0..n {
458            let v_ind = horseshoe_influence(&panels[j], panels[i].control_point, far_field);
459            matrix[i][j] = v_ind.dot(panels[i].normal);
460        }
461        rhs[i] = -v_free.dot(panels[i].normal);
462    }
463
464    // Solve for circulation strengths
465    let (lu, pivot) = hisab::num::lu_decompose(&matrix)
466        .map_err(|e| PavanError::ComputationError(format!("VLM LU decomposition failed: {e}")))?;
467    let gamma = hisab::num::lu_solve(&lu, &pivot, &rhs)
468        .map_err(|e| PavanError::ComputationError(format!("VLM LU solve failed: {e}")))?;
469
470    Ok(vlm_post_process(panels, gamma, wing, v_inf))
471}
472
473/// Solve VLM for multiple angles of attack efficiently.
474///
475/// The influence matrix is LU-decomposed once and reused for each angle.
476#[allow(clippy::needless_range_loop)]
477pub fn solve_multi(
478    panels: &[VlmPanel],
479    wing: &WingGeometry,
480    alphas: &[f64],
481    v_inf: f64,
482) -> Result<Vec<VlmSolution>> {
483    let n = panels.len();
484    if n < MIN_VLM_PANELS {
485        return Err(PavanError::InvalidGeometry(format!(
486            "need at least {MIN_VLM_PANELS} VLM panels, got {n}"
487        )));
488    }
489    if v_inf <= 0.0 {
490        return Err(PavanError::InvalidVelocity(
491            "freestream velocity must be positive".into(),
492        ));
493    }
494
495    let far_field = wing.root_chord * TRAILING_LEG_FACTOR;
496
497    // Build influence matrix (angle-independent)
498    let mut matrix: Vec<Vec<f64>> = vec![vec![0.0; n]; n];
499    for i in 0..n {
500        for j in 0..n {
501            let v_ind = horseshoe_influence(&panels[j], panels[i].control_point, far_field);
502            matrix[i][j] = v_ind.dot(panels[i].normal);
503        }
504    }
505
506    let (lu, pivot) = hisab::num::lu_decompose(&matrix)
507        .map_err(|e| PavanError::ComputationError(format!("VLM LU decomposition failed: {e}")))?;
508
509    let mut results = Vec::with_capacity(alphas.len());
510
511    for &alpha in alphas {
512        let v_free = DVec3::new(v_inf * alpha.cos(), 0.0, v_inf * alpha.sin());
513
514        let rhs: Vec<f64> = panels.iter().map(|p| -v_free.dot(p.normal)).collect();
515
516        let gamma = hisab::num::lu_solve(&lu, &pivot, &rhs)
517            .map_err(|e| PavanError::ComputationError(format!("VLM LU solve failed: {e}")))?;
518
519        results.push(vlm_post_process(panels, gamma, wing, v_inf));
520    }
521
522    Ok(results)
523}
524
525/// Compute Oswald span efficiency from VLM results.
526///
527/// e = CL² / (π · AR · CDi). Returns 1.0 for elliptic loading.
528#[must_use]
529#[inline]
530pub fn oswald_efficiency(cl: f64, cdi: f64, aspect_ratio: f64) -> f64 {
531    if cdi.abs() < f64::EPSILON || aspect_ratio <= 0.0 {
532        return 1.0;
533    }
534    cl * cl / (PI * aspect_ratio * cdi)
535}
536
537#[cfg(test)]
538mod tests {
539    use super::*;
540
541    fn rect_wing(n_span: usize, n_chord: usize) -> (WingGeometry, Vec<VlmPanel>) {
542        // AR=6 rectangular wing: span=6m, chord=1m
543        let wing = WingGeometry::rectangular(6.0, 1.0, n_span, n_chord);
544        let panels = generate_panels(&wing);
545        (wing, panels)
546    }
547
548    // --- Panel generation ---
549
550    #[test]
551    fn panel_count_correct() {
552        let wing = WingGeometry::rectangular(6.0, 1.0, 10, 4);
553        let panels = generate_panels(&wing);
554        assert_eq!(panels.len(), 2 * 10 * 4);
555    }
556
557    #[test]
558    fn panels_have_positive_area() {
559        let (_, panels) = rect_wing(10, 2);
560        for (i, p) in panels.iter().enumerate() {
561            assert!(p.area > 0.0, "panel {i} has non-positive area {}", p.area);
562        }
563    }
564
565    #[test]
566    fn panel_normals_are_unit() {
567        let (_, panels) = rect_wing(10, 2);
568        for (i, p) in panels.iter().enumerate() {
569            let mag = p.normal.length();
570            assert!(
571                (mag - 1.0).abs() < 1e-10,
572                "panel {i} normal not unit: {mag}"
573            );
574        }
575    }
576
577    #[test]
578    fn panel_normals_point_upward() {
579        // For a flat wing in the XY plane, normals should point +Z
580        let (_, panels) = rect_wing(10, 2);
581        for (i, p) in panels.iter().enumerate() {
582            assert!(
583                p.normal.z > 0.5,
584                "panel {i} normal should point up, got z={}",
585                p.normal.z
586            );
587        }
588    }
589
590    #[test]
591    fn control_points_within_wing() {
592        let (wing, panels) = rect_wing(10, 2);
593        let half_span = wing.span * 0.5;
594        for p in &panels {
595            assert!(p.control_point.y.abs() <= half_span + 0.01);
596            assert!(p.control_point.x >= -0.01);
597            assert!(p.control_point.x <= wing.root_chord + 0.01);
598        }
599    }
600
601    // --- Biot-Savart ---
602
603    #[test]
604    fn biot_savart_known_geometry() {
605        // Point at (0.5,1,0) from segment along x-axis (0,0,0)→(1,0,0)
606        // Right-hand rule: dl(+x) × r(+y) = +z
607        let v = biot_savart_segment(
608            DVec3::new(0.5, 1.0, 0.0),
609            DVec3::new(0.0, 0.0, 0.0),
610            DVec3::new(1.0, 0.0, 0.0),
611        );
612        assert!(
613            v.z > 0.0,
614            "Biot-Savart should induce +z velocity, got {v:?}"
615        );
616        assert!(v.x.abs() < 1e-10, "x component should be ~0");
617        assert!(v.y.abs() < 1e-10, "y component should be ~0");
618    }
619
620    #[test]
621    fn biot_savart_point_on_line_returns_zero() {
622        // Point on the vortex line itself → core cutoff returns zero
623        let v = biot_savart_segment(
624            DVec3::new(0.5, 0.0, 0.0),
625            DVec3::new(0.0, 0.0, 0.0),
626            DVec3::new(1.0, 0.0, 0.0),
627        );
628        assert_eq!(v, DVec3::ZERO);
629    }
630
631    // --- Solve correctness ---
632
633    #[test]
634    fn rectangular_wing_ar6_at_5deg() {
635        let (wing, panels) = rect_wing(10, 2);
636        let sol = solve(&panels, &wing, 5.0_f64.to_radians(), 1.0).expect("solve");
637        // Finite wing CL < 2D Cl (0.55). For AR=6, CL ≈ 0.35-0.50
638        assert!(
639            sol.cl > 0.2 && sol.cl < 0.6,
640            "AR=6 wing at 5° should have CL in [0.2, 0.6], got {}",
641            sol.cl
642        );
643    }
644
645    #[test]
646    fn positive_induced_drag() {
647        let (wing, panels) = rect_wing(10, 2);
648        let sol = solve(&panels, &wing, 5.0_f64.to_radians(), 1.0).expect("solve");
649        assert!(
650            sol.cdi > 0.0,
651            "induced drag should be positive, got {}",
652            sol.cdi
653        );
654    }
655
656    #[test]
657    fn symmetric_wing_zero_aoa_zero_lift() {
658        let (wing, panels) = rect_wing(8, 2);
659        let sol = solve(&panels, &wing, 0.0, 1.0).expect("solve");
660        assert!(
661            sol.cl.abs() < 0.01,
662            "symmetric wing at 0° should have CL ≈ 0, got {}",
663            sol.cl
664        );
665    }
666
667    #[test]
668    fn lift_increases_with_aoa() {
669        let (wing, panels) = rect_wing(8, 2);
670        let sol2 = solve(&panels, &wing, 2.0_f64.to_radians(), 1.0).expect("solve");
671        let sol5 = solve(&panels, &wing, 5.0_f64.to_radians(), 1.0).expect("solve");
672        let sol8 = solve(&panels, &wing, 8.0_f64.to_radians(), 1.0).expect("solve");
673        assert!(sol5.cl > sol2.cl);
674        assert!(sol8.cl > sol5.cl);
675    }
676
677    #[test]
678    fn oswald_efficiency_reasonable() {
679        let (wing, panels) = rect_wing(10, 2);
680        let sol = solve(&panels, &wing, 5.0_f64.to_radians(), 1.0).expect("solve");
681        // Rectangular wing: e ≈ 0.7-1.0
682        assert!(
683            sol.oswald_efficiency > 0.5 && sol.oswald_efficiency < 1.5,
684            "Oswald e should be reasonable, got {}",
685            sol.oswald_efficiency
686        );
687    }
688
689    #[test]
690    fn span_loading_has_correct_length() {
691        let (wing, panels) = rect_wing(10, 2);
692        let sol = solve(&panels, &wing, 5.0_f64.to_radians(), 1.0).expect("solve");
693        assert_eq!(sol.span_cl.len(), 2 * 10);
694    }
695
696    #[test]
697    fn tapered_wing_positive_lift() {
698        let wing = WingGeometry::tapered(6.0, 1.5, 0.5, 10, 2);
699        let panels = generate_panels(&wing);
700        let sol = solve(&panels, &wing, 5.0_f64.to_radians(), 1.0).expect("solve");
701        assert!(sol.cl > 0.1, "tapered wing should produce lift");
702    }
703
704    #[test]
705    fn convergence_with_panel_count() {
706        let wing_lo = WingGeometry::rectangular(6.0, 1.0, 5, 1);
707        let wing_hi = WingGeometry::rectangular(6.0, 1.0, 20, 2);
708        let panels_lo = generate_panels(&wing_lo);
709        let panels_hi = generate_panels(&wing_hi);
710
711        let sol_lo = solve(&panels_lo, &wing_lo, 5.0_f64.to_radians(), 1.0).expect("solve");
712        let sol_hi = solve(&panels_hi, &wing_hi, 5.0_f64.to_radians(), 1.0).expect("solve");
713
714        // Both should give roughly similar CL
715        let diff = (sol_lo.cl - sol_hi.cl).abs();
716        assert!(
717            diff < 0.15,
718            "CL should converge: lo={}, hi={}, diff={}",
719            sol_lo.cl,
720            sol_hi.cl,
721            diff
722        );
723    }
724
725    // --- solve_multi ---
726
727    #[test]
728    fn solve_multi_matches_individual() {
729        let (wing, panels) = rect_wing(6, 2);
730        let alphas = [0.0, 3.0_f64.to_radians(), 5.0_f64.to_radians()];
731        let multi = solve_multi(&panels, &wing, &alphas, 1.0).expect("solve_multi");
732
733        for (i, &alpha) in alphas.iter().enumerate() {
734            let single = solve(&panels, &wing, alpha, 1.0).expect("solve");
735            assert!(
736                (multi[i].cl - single.cl).abs() < 1e-10,
737                "CL mismatch at alpha[{i}]"
738            );
739            assert!(
740                (multi[i].cdi - single.cdi).abs() < 1e-10,
741                "CDi mismatch at alpha[{i}]"
742            );
743            assert!(
744                (multi[i].cm - single.cm).abs() < 1e-10,
745                "Cm mismatch at alpha[{i}]"
746            );
747        }
748    }
749
750    #[test]
751    fn solve_multi_empty_alphas() {
752        let (wing, panels) = rect_wing(6, 2);
753        let results = solve_multi(&panels, &wing, &[], 1.0).expect("solve_multi");
754        assert!(results.is_empty());
755    }
756
757    // --- Error cases ---
758
759    #[test]
760    fn too_few_panels_errors() {
761        let wing = WingGeometry::rectangular(6.0, 1.0, 1, 1);
762        let panels = vec![]; // empty
763        assert!(solve(&panels, &wing, 0.0, 1.0).is_err());
764    }
765
766    #[test]
767    fn zero_velocity_errors() {
768        let (wing, panels) = rect_wing(6, 2);
769        assert!(solve(&panels, &wing, 0.0, 0.0).is_err());
770    }
771
772    // --- Geometry presets ---
773
774    #[test]
775    fn rectangular_reference_area() {
776        let wing = WingGeometry::rectangular(6.0, 1.0, 10, 2);
777        assert!((wing.reference_area() - 6.0).abs() < f64::EPSILON);
778    }
779
780    #[test]
781    fn rectangular_aspect_ratio() {
782        let wing = WingGeometry::rectangular(6.0, 1.0, 10, 2);
783        assert!((wing.aspect_ratio() - 6.0).abs() < f64::EPSILON);
784    }
785
786    #[test]
787    fn tapered_reference_area() {
788        let wing = WingGeometry::tapered(10.0, 2.0, 1.0, 10, 2);
789        // S = 10 × (2 + 1) / 2 = 15
790        assert!((wing.reference_area() - 15.0).abs() < f64::EPSILON);
791    }
792
793    #[test]
794    fn oswald_efficiency_fn() {
795        // CL=0.5, CDi=0.005, AR=8 → e = 0.25/(π·8·0.005) = 0.25/0.1257 ≈ 1.99
796        // That's high — just testing the formula works
797        let e = oswald_efficiency(0.5, 0.0133, 6.0);
798        assert!(e > 0.5 && e < 1.5, "e should be reasonable, got {e}");
799    }
800
801    // --- Guard clause coverage ---
802
803    #[test]
804    fn total_panels_count() {
805        let wing = WingGeometry::rectangular(6.0, 1.0, 10, 3);
806        assert_eq!(wing.total_panels(), 60);
807    }
808
809    #[test]
810    fn zero_span_area_and_ar() {
811        let wing = WingGeometry::rectangular(0.0, 1.0, 5, 2);
812        assert_eq!(wing.reference_area(), 0.0);
813        assert_eq!(wing.aspect_ratio(), 0.0);
814    }
815
816    #[test]
817    fn solve_zero_velocity_errors() {
818        let (wing, panels) = rect_wing(6, 2);
819        assert!(solve(&panels, &wing, 0.0, 0.0).is_err());
820    }
821
822    #[test]
823    fn solve_multi_zero_velocity_errors() {
824        let (wing, panels) = rect_wing(6, 2);
825        assert!(solve_multi(&panels, &wing, &[0.0], 0.0).is_err());
826    }
827
828    #[test]
829    fn solve_multi_too_few_panels_errors() {
830        let wing = WingGeometry::rectangular(6.0, 1.0, 1, 1);
831        assert!(solve_multi(&[], &wing, &[0.0], 1.0).is_err());
832    }
833
834    #[test]
835    fn biot_savart_coincident_endpoints() {
836        // a == b → zero-length segment
837        let v = biot_savart_segment(
838            DVec3::new(1.0, 1.0, 0.0),
839            DVec3::new(0.0, 0.0, 0.0),
840            DVec3::new(0.0, 0.0, 0.0),
841        );
842        assert_eq!(v, DVec3::ZERO);
843    }
844
845    #[test]
846    fn biot_savart_point_near_endpoint() {
847        let v = biot_savart_segment(
848            DVec3::new(1e-15, 0.0, 0.0),
849            DVec3::new(0.0, 0.0, 0.0),
850            DVec3::new(1.0, 0.0, 0.0),
851        );
852        // Should return zero due to core cutoff
853        assert_eq!(v, DVec3::ZERO);
854    }
855
856    #[test]
857    fn wing_geometry_serde() {
858        let wing = WingGeometry::rectangular(6.0, 1.0, 10, 2);
859        let json = serde_json::to_string(&wing).expect("serialize");
860        let back: WingGeometry = serde_json::from_str(&json).expect("deserialize");
861        assert!((back.span - wing.span).abs() < f64::EPSILON);
862    }
863
864    #[test]
865    fn vlm_solution_serde() {
866        let (wing, panels) = rect_wing(6, 2);
867        let sol = solve(&panels, &wing, 5.0_f64.to_radians(), 1.0).expect("solve");
868        let json = serde_json::to_string(&sol).expect("serialize");
869        let back: VlmSolution = serde_json::from_str(&json).expect("deserialize");
870        assert!((back.cl - sol.cl).abs() < f64::EPSILON);
871    }
872}