Skip to main content

oxiphysics_geometry/
bspline.rs

1#![allow(clippy::needless_range_loop)]
2// Copyright 2026 COOLJAPAN OU (Team KitaSan)
3// SPDX-License-Identifier: Apache-2.0
4
5//! B-spline and NURBS curves and surfaces.
6//!
7//! Provides:
8//!
9//! - [`BsplineBasis`]: B-spline basis functions N_{i,p}(t) via Cox-de Boor recursion
10//! - [`BsplineCurve`]: Evaluate curve, derivative, arc length, curvature, torsion
11//! - [`BsplineSurface`]: Tensor-product surface, normal, Gaussian/mean curvature
12//! - [`NurbsCurve`]: Rational B-spline, weight update, circle/arc as exact NURBS
13//! - [`NurbsSurface`]: Rational surface, trim curves, surface fitting to point cloud
14//! - [`BsplineFitting`]: Least-squares fitting, chord-length parameterization, knot selection
15
16#![allow(dead_code)]
17#![allow(clippy::too_many_arguments)]
18
19// ============================================================================
20// BsplineBasis
21// ============================================================================
22
23/// B-spline knot vector with associated utilities.
24///
25/// A knot vector is a non-decreasing sequence of real numbers
26/// T = \[t_0, t_1, …, t_{n+p+1}\] where n+1 is the number of basis functions
27/// and p is the polynomial degree.
28#[derive(Debug, Clone)]
29pub struct KnotVector {
30    /// Knot values (non-decreasing).
31    pub knots: Vec<f64>,
32}
33
34impl KnotVector {
35    /// Construct a knot vector from a slice.
36    ///
37    /// # Panics
38    ///
39    /// Panics if the knots are not non-decreasing.
40    pub fn new(knots: Vec<f64>) -> Self {
41        for i in 1..knots.len() {
42            assert!(
43                knots[i] >= knots[i - 1],
44                "knot vector must be non-decreasing"
45            );
46        }
47        Self { knots }
48    }
49
50    /// Construct a uniform (clamped) knot vector for n+1 control points and degree p.
51    ///
52    /// The first and last knots are repeated p+1 times (clamped), and the interior
53    /// knots are uniformly spaced.
54    pub fn clamped_uniform(n_ctrl: usize, degree: usize) -> Self {
55        let n = n_ctrl - 1;
56        let p = degree;
57        let m = n + p + 1;
58        let mut knots = vec![0.0_f64; m + 1];
59        // First p+1 knots = 0
60        for i in 0..=p {
61            knots[i] = 0.0;
62        }
63        // Last p+1 knots = 1
64        for i in (m - p)..=m {
65            knots[i] = 1.0;
66        }
67        // Interior knots
68        let n_interior = m - 2 * p - 1;
69        for j in 1..=n_interior {
70            knots[p + j] = j as f64 / (n_interior + 1) as f64;
71        }
72        Self { knots }
73    }
74
75    /// Find the knot span index i such that t ∈ \[t_i, t_{i+1}).
76    ///
77    /// Uses binary search. For t at the end of the domain, returns the last
78    /// non-zero-length span.
79    pub fn find_span(&self, t: f64, degree: usize, n_ctrl: usize) -> usize {
80        let n = n_ctrl - 1;
81        let p = degree;
82        let knots = &self.knots;
83        // Special case: t == last knot
84        if t >= knots[n + 1] {
85            // Find rightmost span with non-zero length
86            let mut i = n;
87            while i > p && (knots[i] - knots[i + 1]).abs() < 1e-14 {
88                i -= 1;
89            }
90            return i;
91        }
92        if t <= knots[p] {
93            return p;
94        }
95        // Binary search
96        let mut lo = p;
97        let mut hi = n + 1;
98        let mut mid = (lo + hi) / 2;
99        while t < knots[mid] || t >= knots[mid + 1] {
100            if t < knots[mid] {
101                hi = mid;
102            } else {
103                lo = mid;
104            }
105            mid = (lo + hi) / 2;
106        }
107        mid
108    }
109
110    /// Number of knots in the vector.
111    pub fn len(&self) -> usize {
112        self.knots.len()
113    }
114
115    /// Whether the knot vector is empty.
116    pub fn is_empty(&self) -> bool {
117        self.knots.is_empty()
118    }
119
120    /// Return the domain \[a, b\] of the knot vector (a = first knot, b = last knot).
121    pub fn domain(&self) -> (f64, f64) {
122        (
123            *self.knots.first().unwrap_or(&0.0),
124            *self.knots.last().unwrap_or(&1.0),
125        )
126    }
127}
128
129/// B-spline basis functions and their derivatives.
130///
131/// Implements the Cox-de Boor recursion:
132///
133/// N_{i,0}(t) = 1 if t_i ≤ t < t_{i+1}, else 0
134///
135/// N_{i,p}(t) = (t − t_i)/(t_{i+p} − t_i) N_{i,p-1}(t) + (t_{i+p+1} − t)/(t_{i+p+1} − t_{i+1}) N_{i+1,p-1}(t)
136#[derive(Debug, Clone)]
137pub struct BsplineBasis {
138    /// Polynomial degree p.
139    pub degree: usize,
140    /// Knot vector.
141    pub knot_vector: KnotVector,
142    /// Number of control points (n+1).
143    pub n_ctrl: usize,
144}
145
146impl BsplineBasis {
147    /// Create a new B-spline basis.
148    pub fn new(degree: usize, knot_vector: KnotVector, n_ctrl: usize) -> Self {
149        Self {
150            degree,
151            knot_vector,
152            n_ctrl,
153        }
154    }
155
156    /// Create a basis with a uniform clamped knot vector.
157    pub fn clamped(degree: usize, n_ctrl: usize) -> Self {
158        let kv = KnotVector::clamped_uniform(n_ctrl, degree);
159        Self {
160            degree,
161            knot_vector: kv,
162            n_ctrl,
163        }
164    }
165
166    /// Evaluate all non-zero basis functions N_{i-p,p}(t), …, N_{i,p}(t) at parameter t.
167    ///
168    /// Returns an array of length p+1. Uses the triangular table algorithm.
169    pub fn eval_nonzero(&self, t: f64) -> (usize, Vec<f64>) {
170        let p = self.degree;
171        let knots = &self.knot_vector.knots;
172        let span = self.knot_vector.find_span(t, p, self.n_ctrl);
173        let mut n = vec![0.0_f64; p + 1];
174        let mut left = vec![0.0_f64; p + 1];
175        let mut right = vec![0.0_f64; p + 1];
176        n[0] = 1.0;
177        for j in 1..=p {
178            left[j] = t - knots[span + 1 - j];
179            right[j] = knots[span + j] - t;
180            let mut saved = 0.0;
181            for r in 0..j {
182                let temp = n[r] / (right[r + 1] + left[j - r]);
183                n[r] = saved + right[r + 1] * temp;
184                saved = left[j - r] * temp;
185            }
186            n[j] = saved;
187        }
188        (span, n)
189    }
190
191    /// Evaluate all basis functions N_{i,p}(t) for i = 0, …, n_ctrl-1.
192    ///
193    /// Returns a vector of length n_ctrl.
194    pub fn eval_all(&self, t: f64) -> Vec<f64> {
195        let (span, nonzero) = self.eval_nonzero(t);
196        let p = self.degree;
197        let mut result = vec![0.0_f64; self.n_ctrl];
198        for k in 0..=p {
199            if span + k >= p && span + k - p < self.n_ctrl {
200                result[span + k - p] = nonzero[k];
201            }
202        }
203        result
204    }
205
206    /// Evaluate basis function N_{i,p}(t).
207    pub fn eval(&self, i: usize, t: f64) -> f64 {
208        let all = self.eval_all(t);
209        if i < all.len() { all[i] } else { 0.0 }
210    }
211
212    /// Evaluate all non-zero basis functions and their derivatives up to order `deriv`.
213    ///
214    /// Returns a 2D array `ders[k][j]` where k is the derivative order (0..=deriv)
215    /// and j is the basis function index (0..=p) relative to the span.
216    pub fn eval_nonzero_derivs(&self, t: f64, deriv: usize) -> (usize, Vec<Vec<f64>>) {
217        let p = self.degree;
218        let knots = &self.knot_vector.knots;
219        let span = self.knot_vector.find_span(t, p, self.n_ctrl);
220        let d = deriv.min(p);
221        // ndu table
222        let mut ndu = vec![vec![0.0_f64; p + 1]; p + 1];
223        let mut left = vec![0.0_f64; p + 1];
224        let mut right = vec![0.0_f64; p + 1];
225        ndu[0][0] = 1.0;
226        for j in 1..=p {
227            left[j] = t - knots[span + 1 - j];
228            right[j] = knots[span + j] - t;
229            let mut saved = 0.0;
230            for r in 0..j {
231                ndu[j][r] = right[r + 1] + left[j - r];
232                let temp = ndu[r][j - 1] / ndu[j][r];
233                ndu[r][j] = saved + right[r + 1] * temp;
234                saved = left[j - r] * temp;
235            }
236            ndu[j][j] = saved;
237        }
238        let mut ders = vec![vec![0.0_f64; p + 1]; d + 1];
239        for j in 0..=p {
240            ders[0][j] = ndu[j][p];
241        }
242        let mut a = vec![vec![0.0_f64; p + 1]; 2];
243        for r in 0..=p {
244            let mut s1 = 0_usize;
245            let mut s2 = 1_usize;
246            a[0][0] = 1.0;
247            for k in 1..=d {
248                let mut nd = 0.0;
249                let rk = r as isize - k as isize;
250                let pk = p as isize - k as isize;
251                if r >= k {
252                    a[s2][0] = a[s1][0] / ndu[pk as usize + 1][rk as usize];
253                    nd = a[s2][0] * ndu[rk as usize][pk as usize];
254                }
255                let j1 = if rk >= -1 {
256                    1_usize
257                } else {
258                    (-(rk + 1)) as usize
259                };
260                let j2 = if (r as isize - 1) <= pk { k - 1 } else { p - r };
261                for j in j1..=j2 {
262                    let idx2 = rk + j as isize;
263                    if idx2 < 0 {
264                        continue;
265                    }
266                    let idx2 = idx2 as usize;
267                    a[s2][j] = (a[s1][j] - a[s1][j.saturating_sub(1)]) / ndu[pk as usize + 1][idx2];
268                    nd += a[s2][j] * ndu[idx2][pk as usize];
269                }
270                if r <= (p - k) {
271                    a[s2][k] = -a[s1][k - 1] / ndu[pk as usize + 1][r];
272                    nd += a[s2][k] * ndu[r][pk as usize];
273                }
274                ders[k][r] = nd;
275                std::mem::swap(&mut s1, &mut s2);
276            }
277        }
278        let mut rfact = p as f64;
279        for k in 1..=d {
280            for j in 0..=p {
281                ders[k][j] *= rfact;
282            }
283            if k < d {
284                rfact *= (p - k) as f64;
285            }
286        }
287        (span, ders)
288    }
289
290    /// Evaluate the k-th derivative of basis function N_{i,p}(t).
291    pub fn eval_deriv(&self, i: usize, t: f64, k: usize) -> f64 {
292        let p = self.degree;
293        let (span, ders) = self.eval_nonzero_derivs(t, k);
294        let local_idx = i as isize - span as isize;
295        if local_idx < 0 || local_idx > p as isize {
296            0.0
297        } else {
298            ders[k][local_idx as usize]
299        }
300    }
301
302    /// Greville abscissae (average knots) for control point parameterization.
303    ///
304    /// ξ_i = (t_{i+1} + … + t_{i+p}) / p
305    pub fn greville_abscissae(&self) -> Vec<f64> {
306        let p = self.degree;
307        let knots = &self.knot_vector.knots;
308        (0..self.n_ctrl)
309            .map(|i| {
310                if p == 0 {
311                    knots[i]
312                } else {
313                    knots[i + 1..i + p + 1].iter().sum::<f64>() / p as f64
314                }
315            })
316            .collect()
317    }
318}
319
320// ============================================================================
321// BsplineCurve
322// ============================================================================
323
324/// A B-spline curve in 3D space.
325///
326/// C(t) = Σ_{i=0}^{n} N_{i,p}(t) * P_i
327///
328/// where P_i are control points and N_{i,p}(t) are B-spline basis functions.
329#[derive(Debug, Clone)]
330pub struct BsplineCurve {
331    /// B-spline basis.
332    pub basis: BsplineBasis,
333    /// Control points \[P_0, …, P_n\], each \[x, y, z\].
334    pub control_points: Vec<[f64; 3]>,
335}
336
337impl BsplineCurve {
338    /// Create a new B-spline curve.
339    pub fn new(degree: usize, knot_vector: KnotVector, control_points: Vec<[f64; 3]>) -> Self {
340        let n_ctrl = control_points.len();
341        let basis = BsplineBasis::new(degree, knot_vector, n_ctrl);
342        Self {
343            basis,
344            control_points,
345        }
346    }
347
348    /// Create a B-spline curve with a clamped uniform knot vector.
349    pub fn clamped(degree: usize, control_points: Vec<[f64; 3]>) -> Self {
350        let n_ctrl = control_points.len();
351        let basis = BsplineBasis::clamped(degree, n_ctrl);
352        Self {
353            basis,
354            control_points,
355        }
356    }
357
358    /// Evaluate the curve at parameter t: returns \[x, y, z\].
359    pub fn eval(&self, t: f64) -> [f64; 3] {
360        let (span, nonzero) = self.basis.eval_nonzero(t);
361        let p = self.basis.degree;
362        let mut point = [0.0_f64; 3];
363        for k in 0..=p {
364            let idx = span + k - p;
365            if idx < self.control_points.len() {
366                let cp = self.control_points[idx];
367                for d in 0..3 {
368                    point[d] += nonzero[k] * cp[d];
369                }
370            }
371        }
372        point
373    }
374
375    /// Evaluate the k-th derivative of the curve at parameter t.
376    ///
377    /// Returns \[dx^k/dt^k, dy^k/dt^k, dz^k/dt^k\].
378    pub fn eval_deriv(&self, t: f64, k: usize) -> [f64; 3] {
379        let p = self.basis.degree;
380        let (span, ders) = self.basis.eval_nonzero_derivs(t, k);
381        let n_ctrl = self.control_points.len();
382        let mut result = [0.0_f64; 3];
383        for j in 0..=p {
384            let idx = span + j - p;
385            if idx < n_ctrl {
386                for d in 0..3 {
387                    result[d] += ders[k][j] * self.control_points[idx][d];
388                }
389            }
390        }
391        result
392    }
393
394    /// Compute arc length from t_start to t_end using Gaussian quadrature (n_gauss points per interval).
395    pub fn arc_length(&self, t_start: f64, t_end: f64, n_intervals: usize) -> f64 {
396        let h = (t_end - t_start) / n_intervals as f64;
397        // 5-point Gauss-Legendre on [-1,1]
398        let gp = [-0.906180, -0.538469, 0.0, 0.538469, 0.906180];
399        let gw = [0.236927, 0.478629, 0.568889, 0.478629, 0.236927];
400        let mut length = 0.0;
401        for i in 0..n_intervals {
402            let a = t_start + i as f64 * h;
403            let b = a + h;
404            let mid = (a + b) * 0.5;
405            let half = (b - a) * 0.5;
406            for (&xi, &wi) in gp.iter().zip(gw.iter()) {
407                let t = mid + half * xi;
408                let dt = self.eval_deriv(t, 1);
409                let speed = (dt[0] * dt[0] + dt[1] * dt[1] + dt[2] * dt[2]).sqrt();
410                length += wi * half * speed;
411            }
412        }
413        length
414    }
415
416    /// Compute the unit tangent vector T(t) = C'(t)/|C'(t)|.
417    pub fn tangent(&self, t: f64) -> [f64; 3] {
418        let d1 = self.eval_deriv(t, 1);
419        let mag = (d1[0] * d1[0] + d1[1] * d1[1] + d1[2] * d1[2]).sqrt();
420        if mag < 1e-14 {
421            [0.0, 0.0, 0.0]
422        } else {
423            [d1[0] / mag, d1[1] / mag, d1[2] / mag]
424        }
425    }
426
427    /// Compute curvature κ(t) = |C' × C''| / |C'|³.
428    pub fn curvature(&self, t: f64) -> f64 {
429        let d1 = self.eval_deriv(t, 1);
430        let d2 = self.eval_deriv(t, 2);
431        let cross = [
432            d1[1] * d2[2] - d1[2] * d2[1],
433            d1[2] * d2[0] - d1[0] * d2[2],
434            d1[0] * d2[1] - d1[1] * d2[0],
435        ];
436        let cross_mag = (cross[0] * cross[0] + cross[1] * cross[1] + cross[2] * cross[2]).sqrt();
437        let d1_mag = (d1[0] * d1[0] + d1[1] * d1[1] + d1[2] * d1[2]).sqrt();
438        if d1_mag < 1e-14 {
439            0.0
440        } else {
441            cross_mag / d1_mag.powi(3)
442        }
443    }
444
445    /// Compute torsion τ(t) = (C' × C'')·C''' / |C' × C''|².
446    pub fn torsion(&self, t: f64) -> f64 {
447        let d1 = self.eval_deriv(t, 1);
448        let d2 = self.eval_deriv(t, 2);
449        let d3 = self.eval_deriv(t, 3);
450        let cross = [
451            d1[1] * d2[2] - d1[2] * d2[1],
452            d1[2] * d2[0] - d1[0] * d2[2],
453            d1[0] * d2[1] - d1[1] * d2[0],
454        ];
455        let cross_mag2 = cross[0] * cross[0] + cross[1] * cross[1] + cross[2] * cross[2];
456        if cross_mag2 < 1e-28 {
457            return 0.0;
458        }
459        let dot_d3 = cross[0] * d3[0] + cross[1] * d3[1] + cross[2] * d3[2];
460        dot_d3 / cross_mag2
461    }
462
463    /// Compute the principal normal vector N(t) = T'(t)/|T'(t)|.
464    pub fn principal_normal(&self, t: f64) -> [f64; 3] {
465        let kappa = self.curvature(t);
466        if kappa < 1e-14 {
467            return [0.0, 0.0, 0.0];
468        }
469        let d1 = self.eval_deriv(t, 1);
470        let d2 = self.eval_deriv(t, 2);
471        let d1_mag = (d1[0] * d1[0] + d1[1] * d1[1] + d1[2] * d1[2]).sqrt();
472        let d1_mag3 = d1_mag.powi(3);
473        // N = (d2 * |d1|^2 - d1 * (d1·d2)) / (kappa * |d1|^3)
474        let d1_dot_d2 = d1[0] * d2[0] + d1[1] * d2[1] + d1[2] * d2[2];
475        let d1_mag2 = d1_mag * d1_mag;
476        let mut normal = [0.0_f64; 3];
477        for i in 0..3 {
478            normal[i] = (d2[i] * d1_mag2 - d1[i] * d1_dot_d2) / (kappa * d1_mag3);
479        }
480        normal
481    }
482
483    /// Evaluate the curve at n_points uniformly spaced parameter values.
484    ///
485    /// Returns a vector of \[x, y, z\] points.
486    pub fn sample(&self, n_points: usize) -> Vec<[f64; 3]> {
487        let (t0, t1) = self.basis.knot_vector.domain();
488        (0..n_points)
489            .map(|i| {
490                let t = t0 + (t1 - t0) * i as f64 / (n_points - 1).max(1) as f64;
491                self.eval(t)
492            })
493            .collect()
494    }
495
496    /// Compute the bounding box of the curve (sampled at n_samples points).
497    ///
498    /// Returns (\[x_min, y_min, z_min\], \[x_max, y_max, z_max\]).
499    pub fn bounding_box(&self, n_samples: usize) -> ([f64; 3], [f64; 3]) {
500        let pts = self.sample(n_samples);
501        let mut lo = pts[0];
502        let mut hi = pts[0];
503        for p in &pts {
504            for d in 0..3 {
505                lo[d] = lo[d].min(p[d]);
506                hi[d] = hi[d].max(p[d]);
507            }
508        }
509        (lo, hi)
510    }
511}
512
513// ============================================================================
514// BsplineSurface
515// ============================================================================
516
517/// A tensor-product B-spline surface in 3D space.
518///
519/// S(u,v) = Σ_i Σ_j N_{i,p}(u) N_{j,q}(v) P_{ij}
520#[derive(Debug, Clone)]
521pub struct BsplineSurface {
522    /// Basis in the u direction.
523    pub basis_u: BsplineBasis,
524    /// Basis in the v direction.
525    pub basis_v: BsplineBasis,
526    /// Control point grid \[i\]\[j\] = \[x, y, z\], row-major (u varies fastest).
527    pub control_points: Vec<Vec<[f64; 3]>>,
528}
529
530impl BsplineSurface {
531    /// Create a new B-spline surface.
532    pub fn new(
533        degree_u: usize,
534        degree_v: usize,
535        knot_u: KnotVector,
536        knot_v: KnotVector,
537        control_points: Vec<Vec<[f64; 3]>>,
538    ) -> Self {
539        let n_u = control_points.len();
540        let n_v = control_points[0].len();
541        Self {
542            basis_u: BsplineBasis::new(degree_u, knot_u, n_u),
543            basis_v: BsplineBasis::new(degree_v, knot_v, n_v),
544            control_points,
545        }
546    }
547
548    /// Create a bilinear patch (p=q=1) from four corner points.
549    pub fn bilinear_patch(p00: [f64; 3], p10: [f64; 3], p01: [f64; 3], p11: [f64; 3]) -> Self {
550        let kv = KnotVector::new(vec![0.0, 0.0, 1.0, 1.0]);
551        let cps = vec![vec![p00, p01], vec![p10, p11]];
552        Self {
553            basis_u: BsplineBasis::new(1, kv.clone(), 2),
554            basis_v: BsplineBasis::new(1, kv, 2),
555            control_points: cps,
556        }
557    }
558
559    /// Evaluate the surface at (u, v): returns \[x, y, z\].
560    pub fn eval(&self, u: f64, v: f64) -> [f64; 3] {
561        let (span_u, nu) = self.basis_u.eval_nonzero(u);
562        let (span_v, nv) = self.basis_v.eval_nonzero(v);
563        let p = self.basis_u.degree;
564        let q = self.basis_v.degree;
565        let n_u = self.control_points.len();
566        let n_v = if n_u > 0 {
567            self.control_points[0].len()
568        } else {
569            0
570        };
571        let mut point = [0.0_f64; 3];
572        for k in 0..=p {
573            let iu = span_u + k - p;
574            if iu >= n_u {
575                continue;
576            }
577            for l in 0..=q {
578                let iv = span_v + l - q;
579                if iv >= n_v {
580                    continue;
581                }
582                let w = nu[k] * nv[l];
583                let cp = self.control_points[iu][iv];
584                for d in 0..3 {
585                    point[d] += w * cp[d];
586                }
587            }
588        }
589        point
590    }
591
592    /// Evaluate partial derivative ∂S/∂u at (u, v).
593    pub fn eval_du(&self, u: f64, v: f64) -> [f64; 3] {
594        let p = self.basis_u.degree;
595        let q = self.basis_v.degree;
596        let (span_u, ders_u) = self.basis_u.eval_nonzero_derivs(u, 1);
597        let (span_v, nv) = self.basis_v.eval_nonzero(v);
598        let n_u = self.control_points.len();
599        let n_v = if n_u > 0 {
600            self.control_points[0].len()
601        } else {
602            0
603        };
604        let mut result = [0.0_f64; 3];
605        for k in 0..=p {
606            let iu = span_u + k - p;
607            if iu >= n_u {
608                continue;
609            }
610            for l in 0..=q {
611                let iv = span_v + l - q;
612                if iv >= n_v {
613                    continue;
614                }
615                let w = ders_u[1][k] * nv[l];
616                let cp = self.control_points[iu][iv];
617                for d in 0..3 {
618                    result[d] += w * cp[d];
619                }
620            }
621        }
622        result
623    }
624
625    /// Evaluate partial derivative ∂S/∂v at (u, v).
626    pub fn eval_dv(&self, u: f64, v: f64) -> [f64; 3] {
627        let p = self.basis_u.degree;
628        let q = self.basis_v.degree;
629        let (span_u, nu) = self.basis_u.eval_nonzero(u);
630        let (span_v, ders_v) = self.basis_v.eval_nonzero_derivs(v, 1);
631        let n_u = self.control_points.len();
632        let n_v = if n_u > 0 {
633            self.control_points[0].len()
634        } else {
635            0
636        };
637        let mut result = [0.0_f64; 3];
638        for k in 0..=p {
639            let iu = span_u + k - p;
640            if iu >= n_u {
641                continue;
642            }
643            for l in 0..=q {
644                let iv = span_v + l - q;
645                if iv >= n_v {
646                    continue;
647                }
648                let w = nu[k] * ders_v[1][l];
649                let cp = self.control_points[iu][iv];
650                for d in 0..3 {
651                    result[d] += w * cp[d];
652                }
653            }
654        }
655        result
656    }
657
658    /// Compute the unit surface normal at (u, v).
659    ///
660    /// N = (∂S/∂u × ∂S/∂v) / |∂S/∂u × ∂S/∂v|
661    pub fn normal(&self, u: f64, v: f64) -> [f64; 3] {
662        let du = self.eval_du(u, v);
663        let dv = self.eval_dv(u, v);
664        let cross = [
665            du[1] * dv[2] - du[2] * dv[1],
666            du[2] * dv[0] - du[0] * dv[2],
667            du[0] * dv[1] - du[1] * dv[0],
668        ];
669        let mag = (cross[0] * cross[0] + cross[1] * cross[1] + cross[2] * cross[2]).sqrt();
670        if mag < 1e-14 {
671            [0.0, 0.0, 1.0]
672        } else {
673            [cross[0] / mag, cross[1] / mag, cross[2] / mag]
674        }
675    }
676
677    /// Compute the first fundamental form coefficients E, F, G at (u, v).
678    ///
679    /// E = ∂S/∂u · ∂S/∂u, F = ∂S/∂u · ∂S/∂v, G = ∂S/∂v · ∂S/∂v
680    pub fn first_fundamental_form(&self, u: f64, v: f64) -> (f64, f64, f64) {
681        let su = self.eval_du(u, v);
682        let sv = self.eval_dv(u, v);
683        let e = su[0] * su[0] + su[1] * su[1] + su[2] * su[2];
684        let f = su[0] * sv[0] + su[1] * sv[1] + su[2] * sv[2];
685        let g = sv[0] * sv[0] + sv[1] * sv[1] + sv[2] * sv[2];
686        (e, f, g)
687    }
688
689    /// Compute Gaussian curvature K = (LN - M²)/(EG - F²) at (u, v).
690    ///
691    /// Uses finite differences for second derivatives.
692    pub fn gaussian_curvature(&self, u: f64, v: f64) -> f64 {
693        let h = 1e-5;
694        let (e, f, g) = self.first_fundamental_form(u, v);
695        let egf2 = e * g - f * f;
696        if egf2.abs() < 1e-20 {
697            return 0.0;
698        }
699        // Second derivatives via finite differences
700        let suu = finite_diff_3d(|uu| self.eval_du(uu, v), u, h);
701        let svv = finite_diff_3d(|vv| self.eval_dv(u, vv), v, h);
702        let suv = finite_diff_3d(|uu| self.eval_dv(uu, v), u, h);
703        let n = self.normal(u, v);
704        let l = dot3(suu, n);
705        let m_coef = dot3(suv, n);
706        let nm = dot3(svv, n);
707        (l * nm - m_coef * m_coef) / egf2
708    }
709
710    /// Compute mean curvature H = (EN - 2FM + GL)/(2(EG - F²)) at (u, v).
711    pub fn mean_curvature(&self, u: f64, v: f64) -> f64 {
712        let h = 1e-5;
713        let (e, f, g) = self.first_fundamental_form(u, v);
714        let egf2 = e * g - f * f;
715        if egf2.abs() < 1e-20 {
716            return 0.0;
717        }
718        let suu = finite_diff_3d(|uu| self.eval_du(uu, v), u, h);
719        let svv = finite_diff_3d(|vv| self.eval_dv(u, vv), v, h);
720        let suv = finite_diff_3d(|uu| self.eval_dv(uu, v), u, h);
721        let n = self.normal(u, v);
722        let l = dot3(suu, n);
723        let m_coef = dot3(suv, n);
724        let nm = dot3(svv, n);
725        (e * nm - 2.0 * f * m_coef + g * l) / (2.0 * egf2)
726    }
727
728    /// Sample the surface at (nu × nv) parameter values.
729    pub fn sample(&self, nu: usize, nv: usize) -> Vec<Vec<[f64; 3]>> {
730        let (u0, u1) = self.basis_u.knot_vector.domain();
731        let (v0, v1) = self.basis_v.knot_vector.domain();
732        (0..nu)
733            .map(|i| {
734                let u = u0 + (u1 - u0) * i as f64 / (nu - 1).max(1) as f64;
735                (0..nv)
736                    .map(|j| {
737                        let v = v0 + (v1 - v0) * j as f64 / (nv - 1).max(1) as f64;
738                        self.eval(u, v)
739                    })
740                    .collect()
741            })
742            .collect()
743    }
744    /// Compute Gaussian and mean curvature at (u, v).
745    ///
746    /// Returns `(K_gauss, H_mean)`.
747    pub fn compute_curvature(&self, u: f64, v: f64) -> (f64, f64) {
748        (self.gaussian_curvature(u, v), self.mean_curvature(u, v))
749    }
750    /// Refine the surface by inserting new knots in u and v directions.
751    ///
752    /// Uses the Oslo algorithm (knot insertion) to preserve geometry.
753    pub fn refine_knots(&self, new_knots_u: &[f64], new_knots_v: &[f64]) -> BsplineSurface {
754        // Phase 1: Insert u knots
755        let n_u = self.control_points.len();
756        let n_v = if n_u > 0 {
757            self.control_points[0].len()
758        } else {
759            0
760        };
761        let p = self.basis_u.degree;
762        let q = self.basis_v.degree;
763        // Refine in u: for each column j, refine the curve defined by
764        // control_points[0..n_u][j] along the u knot vector
765        let mut refined_u_knots = self.basis_u.knot_vector.knots.clone();
766        for &k in new_knots_u {
767            refined_u_knots.push(k);
768        }
769        refined_u_knots.sort_by(|a, b| a.partial_cmp(b).unwrap_or(std::cmp::Ordering::Equal));
770        // Recompute control points after u refinement via Boehm's algorithm
771        let mut new_cp = self.control_points.clone();
772        for &knot in new_knots_u {
773            let n_cur = new_cp.len();
774            let orig_u_knots = &self.basis_u.knot_vector.knots;
775            // Find span in the ORIGINAL knot vector
776            let span = orig_u_knots
777                .iter()
778                .position(|&k| k > knot)
779                .map(|pos| pos.saturating_sub(1))
780                .unwrap_or(n_cur.saturating_sub(1))
781                .min(n_cur);
782            let mut inserted = Vec::with_capacity(n_cur + 1);
783            for i in 0..n_cur + 1 {
784                let mut row = Vec::with_capacity(n_v);
785                for j in 0..n_v {
786                    if i < span.saturating_sub(p) + 1 {
787                        row.push(new_cp[i][j]);
788                    } else if i > span {
789                        row.push(new_cp[i - 1][j]);
790                    } else {
791                        let old_knots = &self.basis_u.knot_vector.knots;
792                        let ki = if i < old_knots.len() {
793                            old_knots[i]
794                        } else {
795                            1.0
796                        };
797                        let ki_p = if i + p < old_knots.len() {
798                            old_knots[i + p]
799                        } else {
800                            1.0
801                        };
802                        let denom = ki_p - ki;
803                        let alpha = if denom.abs() > 1e-14 {
804                            (knot - ki) / denom
805                        } else {
806                            0.5
807                        };
808                        let p0 = if i > 0 {
809                            new_cp[i - 1][j]
810                        } else {
811                            new_cp[0][j]
812                        };
813                        let p1 = if i < n_cur {
814                            new_cp[i][j]
815                        } else {
816                            new_cp[n_cur - 1][j]
817                        };
818                        row.push([
819                            (1.0 - alpha) * p0[0] + alpha * p1[0],
820                            (1.0 - alpha) * p0[1] + alpha * p1[1],
821                            (1.0 - alpha) * p0[2] + alpha * p1[2],
822                        ]);
823                    }
824                }
825                inserted.push(row);
826            }
827            new_cp = inserted;
828        }
829        // Phase 2: Insert v knots using Boehm's algorithm (mirroring u-direction)
830        let n_u_new = new_cp.len();
831        let mut new_cp_v = new_cp;
832        let mut refined_v_knots = self.basis_v.knot_vector.knots.clone();
833        for &k in new_knots_v {
834            refined_v_knots.push(k);
835        }
836        refined_v_knots.sort_by(|a, b| a.partial_cmp(b).unwrap_or(std::cmp::Ordering::Equal));
837        for &knot in new_knots_v {
838            let n_v_cur = if !new_cp_v.is_empty() {
839                new_cp_v[0].len()
840            } else {
841                0
842            };
843            let orig_v_knots = &self.basis_v.knot_vector.knots;
844            // Find span in the ORIGINAL v knot vector
845            let span = orig_v_knots
846                .iter()
847                .position(|&k| k > knot)
848                .map(|pos| pos.saturating_sub(1))
849                .unwrap_or(n_v_cur.saturating_sub(1))
850                .min(n_v_cur);
851            let mut new_rows = Vec::with_capacity(n_u_new);
852            for i in 0..n_u_new {
853                let old_row = &new_cp_v[i];
854                let n_cur = old_row.len();
855                let mut new_row = Vec::with_capacity(n_cur + 1);
856                for j in 0..n_cur + 1 {
857                    if j < span.saturating_sub(q) + 1 {
858                        new_row.push(old_row[j]);
859                    } else if j > span {
860                        new_row.push(old_row[j - 1]);
861                    } else {
862                        let old_v_knots = &self.basis_v.knot_vector.knots;
863                        let kj = if j < old_v_knots.len() {
864                            old_v_knots[j]
865                        } else {
866                            1.0
867                        };
868                        let kj_q = if j + q < old_v_knots.len() {
869                            old_v_knots[j + q]
870                        } else {
871                            1.0
872                        };
873                        let denom = kj_q - kj;
874                        let alpha = if denom.abs() > 1e-14 {
875                            (knot - kj) / denom
876                        } else {
877                            0.5
878                        };
879                        let p0 = if j > 0 { old_row[j - 1] } else { old_row[0] };
880                        let p1 = if j < n_cur {
881                            old_row[j]
882                        } else {
883                            old_row[n_cur - 1]
884                        };
885                        new_row.push([
886                            (1.0 - alpha) * p0[0] + alpha * p1[0],
887                            (1.0 - alpha) * p0[1] + alpha * p1[1],
888                            (1.0 - alpha) * p0[2] + alpha * p1[2],
889                        ]);
890                    }
891                }
892                new_rows.push(new_row);
893            }
894            new_cp_v = new_rows;
895        }
896        let final_n_u = new_cp_v.len();
897        let final_n_v = if final_n_u > 0 { new_cp_v[0].len() } else { 0 };
898        // Build new knot vectors
899        let mut ku = self.basis_u.knot_vector.knots.clone();
900        for &k in new_knots_u {
901            ku.push(k);
902        }
903        ku.sort_by(|a, b| a.partial_cmp(b).unwrap_or(std::cmp::Ordering::Equal));
904        let mut kv = self.basis_v.knot_vector.knots.clone();
905        for &k in new_knots_v {
906            kv.push(k);
907        }
908        kv.sort_by(|a, b| a.partial_cmp(b).unwrap_or(std::cmp::Ordering::Equal));
909        // Ensure knot vectors have correct length (n_ctrl + degree + 1)
910        while ku.len() < final_n_u + p + 1 {
911            ku.push(1.0);
912        }
913        while kv.len() < final_n_v + q + 1 {
914            kv.push(1.0);
915        }
916        BsplineSurface::new(p, q, KnotVector::new(ku), KnotVector::new(kv), new_cp_v)
917    }
918}
919
920// ============================================================================
921// NurbsCurve
922// ============================================================================
923
924/// A rational B-spline (NURBS) curve in 3D space.
925///
926/// C(t) = Σ w_i N_{i,p}(t) P_i / Σ w_i N_{i,p}(t)
927#[derive(Debug, Clone)]
928pub struct NurbsCurve {
929    /// Underlying B-spline basis.
930    pub basis: BsplineBasis,
931    /// Control points \[P_0, …, P_n\], each \[x, y, z\].
932    pub control_points: Vec<[f64; 3]>,
933    /// Weights w_i > 0.
934    pub weights: Vec<f64>,
935}
936
937impl NurbsCurve {
938    /// Create a new NURBS curve.
939    pub fn new(
940        degree: usize,
941        knot_vector: KnotVector,
942        control_points: Vec<[f64; 3]>,
943        weights: Vec<f64>,
944    ) -> Self {
945        let n_ctrl = control_points.len();
946        assert_eq!(
947            weights.len(),
948            n_ctrl,
949            "weights and control points must have the same length"
950        );
951        let basis = BsplineBasis::new(degree, knot_vector, n_ctrl);
952        Self {
953            basis,
954            control_points,
955            weights,
956        }
957    }
958
959    /// Create a NURBS circle arc from center, radius, start and end angles (radians).
960    ///
961    /// Uses the standard 9-point representation for a full circle (degree 2).
962    /// For arcs smaller than 2π, a simplified quadratic NURBS is returned.
963    pub fn circle_arc(center: [f64; 3], radius: f64, start_angle: f64, end_angle: f64) -> Self {
964        // Use a simple quadratic NURBS arc (3 control points for ≤ 180°, etc.)
965        // For generality, build a degree-2 NURBS with the standard construction.
966        let sweep = end_angle - start_angle;
967        let n_arcs = if sweep.abs() <= std::f64::consts::FRAC_PI_2 {
968            1
969        } else if sweep.abs() <= std::f64::consts::PI {
970            2
971        } else if sweep.abs() <= 3.0 * std::f64::consts::FRAC_PI_2 {
972            3
973        } else {
974            4
975        };
976        let n_pts = 2 * n_arcs + 1;
977        let mut ctrl_pts = Vec::with_capacity(n_pts);
978        let mut weights = Vec::with_capacity(n_pts);
979        let d_theta = sweep / n_arcs as f64;
980        let w1 = (d_theta / 2.0).cos();
981        let mut angle = start_angle;
982        // First point on arc
983        ctrl_pts.push([
984            center[0] + radius * angle.cos(),
985            center[1] + radius * angle.sin(),
986            center[2],
987        ]);
988        weights.push(1.0);
989        for _i in 0..n_arcs {
990            let a_mid = angle + d_theta * 0.5;
991            let a_end = angle + d_theta;
992            // Mid control point (off-arc)
993            ctrl_pts.push([
994                center[0] + radius / w1 * a_mid.cos(),
995                center[1] + radius / w1 * a_mid.sin(),
996                center[2],
997            ]);
998            weights.push(w1);
999            // End point on arc
1000            ctrl_pts.push([
1001                center[0] + radius * a_end.cos(),
1002                center[1] + radius * a_end.sin(),
1003                center[2],
1004            ]);
1005            weights.push(1.0);
1006            angle = a_end;
1007        }
1008        // Build clamped knot vector for degree 2
1009        let mut knots = vec![0.0_f64; 3]; // first 3
1010        let knot_step = 1.0 / n_arcs as f64;
1011        for i in 1..n_arcs {
1012            knots.push(i as f64 * knot_step);
1013            knots.push(i as f64 * knot_step);
1014        }
1015        knots.extend(std::iter::repeat_n(1.0_f64, 3));
1016        let kv = KnotVector::new(knots);
1017        Self::new(2, kv, ctrl_pts, weights)
1018    }
1019
1020    /// Evaluate the NURBS curve at parameter t: returns \[x, y, z\].
1021    pub fn eval(&self, t: f64) -> [f64; 3] {
1022        let (span, nonzero) = self.basis.eval_nonzero(t);
1023        let p = self.basis.degree;
1024        let n_ctrl = self.control_points.len();
1025        let mut num = [0.0_f64; 3];
1026        let mut denom = 0.0_f64;
1027        for k in 0..=p {
1028            let idx = span + k - p;
1029            if idx >= n_ctrl {
1030                continue;
1031            }
1032            let wn = self.weights[idx] * nonzero[k];
1033            denom += wn;
1034            let cp = self.control_points[idx];
1035            for d in 0..3 {
1036                num[d] += wn * cp[d];
1037            }
1038        }
1039        if denom.abs() < 1e-30 {
1040            return self.control_points[0];
1041        }
1042        [num[0] / denom, num[1] / denom, num[2] / denom]
1043    }
1044
1045    /// Evaluate the first derivative of the NURBS curve at t.
1046    ///
1047    /// Uses the quotient rule: C'(t) = (W'(t) P̃(t) - W(t)^2 P̃'(t)) / W(t)^2
1048    pub fn eval_deriv(&self, t: f64) -> [f64; 3] {
1049        let h = 1e-7;
1050        let t0 = (t - h).max(self.basis.knot_vector.knots[0]);
1051        let t1 = (t + h).min(*self.basis.knot_vector.knots.last().unwrap_or(&1.0));
1052        let p0 = self.eval(t0);
1053        let p1 = self.eval(t1);
1054        let dt = t1 - t0;
1055        [
1056            (p1[0] - p0[0]) / dt,
1057            (p1[1] - p0[1]) / dt,
1058            (p1[2] - p0[2]) / dt,
1059        ]
1060    }
1061
1062    /// Update a control point and its weight.
1063    pub fn update_control_point(&mut self, i: usize, point: [f64; 3], weight: f64) {
1064        self.control_points[i] = point;
1065        self.weights[i] = weight;
1066    }
1067
1068    /// Compute arc length using Gaussian quadrature.
1069    pub fn arc_length(&self, t_start: f64, t_end: f64, n_intervals: usize) -> f64 {
1070        let h = (t_end - t_start) / n_intervals as f64;
1071        let gp = [-0.906180, -0.538469, 0.0, 0.538469, 0.906180];
1072        let gw = [0.236927, 0.478629, 0.568889, 0.478629, 0.236927];
1073        let mut length = 0.0;
1074        for i in 0..n_intervals {
1075            let a = t_start + i as f64 * h;
1076            let b = a + h;
1077            let mid = (a + b) * 0.5;
1078            let half = (b - a) * 0.5;
1079            for (&xi, &wi) in gp.iter().zip(gw.iter()) {
1080                let t = mid + half * xi;
1081                let dt = self.eval_deriv(t);
1082                let speed = (dt[0] * dt[0] + dt[1] * dt[1] + dt[2] * dt[2]).sqrt();
1083                length += wi * half * speed;
1084            }
1085        }
1086        length
1087    }
1088
1089    /// Compute curvature using the B-spline curve's curvature formula (numerical).
1090    pub fn curvature(&self, t: f64) -> f64 {
1091        let h = 1e-5;
1092        let t_min = self.basis.knot_vector.knots[0];
1093        let t_max = *self.basis.knot_vector.knots.last().unwrap_or(&1.0);
1094        let t0 = (t - h).max(t_min);
1095        let t1 = (t + h).min(t_max);
1096        let dt = t1 - t0;
1097        if dt < 1e-14 {
1098            return 0.0;
1099        }
1100        let d1 = self.eval_deriv(t);
1101        let p0 = self.eval(t0);
1102        let p1 = self.eval(t);
1103        let p2 = self.eval(t1);
1104        let d2 = [
1105            (p2[0] - 2.0 * p1[0] + p0[0]) / (h * h),
1106            (p2[1] - 2.0 * p1[1] + p0[1]) / (h * h),
1107            (p2[2] - 2.0 * p1[2] + p0[2]) / (h * h),
1108        ];
1109        let cross = [
1110            d1[1] * d2[2] - d1[2] * d2[1],
1111            d1[2] * d2[0] - d1[0] * d2[2],
1112            d1[0] * d2[1] - d1[1] * d2[0],
1113        ];
1114        let cross_mag = (cross[0] * cross[0] + cross[1] * cross[1] + cross[2] * cross[2]).sqrt();
1115        let d1_mag = (d1[0] * d1[0] + d1[1] * d1[1] + d1[2] * d1[2]).sqrt();
1116        if d1_mag < 1e-14 {
1117            0.0
1118        } else {
1119            cross_mag / d1_mag.powi(3)
1120        }
1121    }
1122
1123    /// Sample the NURBS curve at n_points uniformly spaced parameter values.
1124    pub fn sample(&self, n_points: usize) -> Vec<[f64; 3]> {
1125        let (t0, t1) = self.basis.knot_vector.domain();
1126        (0..n_points)
1127            .map(|i| {
1128                let t = t0 + (t1 - t0) * i as f64 / (n_points - 1).max(1) as f64;
1129                self.eval(t)
1130            })
1131            .collect()
1132    }
1133    /// Evaluate the B-spline basis function N_{i,p}(t) (Cox-de Boor recursion).
1134    pub fn b_spline_basis(i: usize, p: usize, t: f64, knots: &[f64]) -> f64 {
1135        if p == 0 {
1136            return if knots[i] <= t && t < knots[i + 1] {
1137                1.0
1138            } else if (t - knots[knots.len() - 1]).abs() < 1e-14
1139                && (knots[i + 1] - knots[knots.len() - 1]).abs() < 1e-14
1140            {
1141                // Handle t == last knot for the last basis function
1142                1.0
1143            } else {
1144                0.0
1145            };
1146        }
1147        let left_denom = knots[i + p] - knots[i];
1148        let left = if left_denom.abs() > 1e-14 {
1149            (t - knots[i]) / left_denom * Self::b_spline_basis(i, p - 1, t, knots)
1150        } else {
1151            0.0
1152        };
1153        let right_denom = knots[i + p + 1] - knots[i + 1];
1154        let right = if right_denom.abs() > 1e-14 {
1155            (knots[i + p + 1] - t) / right_denom * Self::b_spline_basis(i + 1, p - 1, t, knots)
1156        } else {
1157            0.0
1158        };
1159        left + right
1160    }
1161    /// Create a NURBS curve from 3D control points, weights, and degree.
1162    ///
1163    /// Automatically creates a clamped uniform knot vector.
1164    pub fn from_points_and_weights(
1165        points: Vec<[f64; 3]>,
1166        weights: Vec<f64>,
1167        degree: usize,
1168    ) -> Self {
1169        let n = points.len();
1170        let m = n + degree + 1;
1171        let mut knots = vec![0.0_f64; degree + 1];
1172        let interior = m - 2 * (degree + 1);
1173        for i in 1..=interior {
1174            knots.push(i as f64 / (interior + 1) as f64);
1175        }
1176        knots.extend(vec![1.0_f64; degree + 1]);
1177        let kv = KnotVector::new(knots);
1178        Self::new(degree, kv, points, weights)
1179    }
1180    /// Create a NURBS circle in the XY plane with the given radius (degree 2, 9 points).
1181    pub fn circle(radius: f64) -> Self {
1182        let w = std::f64::consts::FRAC_1_SQRT_2; // cos(π/4) ≈ 0.7071
1183        let r = radius;
1184        let ctrl = vec![
1185            [r, 0.0, 0.0],
1186            [r, r, 0.0],
1187            [0.0, r, 0.0],
1188            [-r, r, 0.0],
1189            [-r, 0.0, 0.0],
1190            [-r, -r, 0.0],
1191            [0.0, -r, 0.0],
1192            [r, -r, 0.0],
1193            [r, 0.0, 0.0],
1194        ];
1195        let weights = vec![1.0, w, 1.0, w, 1.0, w, 1.0, w, 1.0];
1196        let knots = KnotVector::new(vec![
1197            0.0, 0.0, 0.0, 0.25, 0.25, 0.5, 0.5, 0.75, 0.75, 1.0, 1.0, 1.0,
1198        ]);
1199        Self::new(2, knots, ctrl, weights)
1200    }
1201    /// Compute arc-length reparametrization: return `n_output` uniformly-spaced
1202    /// (by arc length) sample points.
1203    pub fn compute_arc_length_reparametrize(
1204        &self,
1205        n_intervals: usize,
1206        n_output: usize,
1207    ) -> Vec<[f64; 3]> {
1208        let (arc_lengths, params) = self.arc_length_table(n_intervals);
1209        let total_len = *arc_lengths.last().unwrap_or(&0.0);
1210        if total_len < 1e-14 || n_output == 0 {
1211            return Vec::new();
1212        }
1213        let mut result = Vec::with_capacity(n_output);
1214        for k in 0..n_output {
1215            let target = total_len * k as f64 / (n_output - 1).max(1) as f64;
1216            // Binary search for the parameter corresponding to this arc length
1217            let idx = arc_lengths
1218                .partition_point(|&s| s < target)
1219                .min(arc_lengths.len() - 1);
1220            let t = if idx == 0 {
1221                params[0]
1222            } else if (arc_lengths[idx] - arc_lengths[idx - 1]).abs() < 1e-14 {
1223                params[idx]
1224            } else {
1225                let frac =
1226                    (target - arc_lengths[idx - 1]) / (arc_lengths[idx] - arc_lengths[idx - 1]);
1227                params[idx - 1] + frac * (params[idx] - params[idx - 1])
1228            };
1229            result.push(self.eval(t));
1230        }
1231        result
1232    }
1233    /// Build an arc-length look-up table with `n` samples.
1234    ///
1235    /// Returns `(cumulative_arc_lengths, corresponding_parameters)`.
1236    pub fn arc_length_table(&self, n: usize) -> (Vec<f64>, Vec<f64>) {
1237        let (t0, t1) = self.basis.knot_vector.domain();
1238        let mut arc_lengths = Vec::with_capacity(n + 1);
1239        let mut params = Vec::with_capacity(n + 1);
1240        let mut cumulative = 0.0_f64;
1241        let mut prev_pt = self.eval(t0);
1242        arc_lengths.push(0.0);
1243        params.push(t0);
1244        for i in 1..=n {
1245            let t = t0 + (t1 - t0) * i as f64 / n as f64;
1246            let pt = self.eval(t);
1247            let dx = pt[0] - prev_pt[0];
1248            let dy = pt[1] - prev_pt[1];
1249            let dz = pt[2] - prev_pt[2];
1250            cumulative += (dx * dx + dy * dy + dz * dz).sqrt();
1251            arc_lengths.push(cumulative);
1252            params.push(t);
1253            prev_pt = pt;
1254        }
1255        (arc_lengths, params)
1256    }
1257}
1258
1259// ============================================================================
1260// NurbsSurface
1261// ============================================================================
1262
1263/// A rational tensor-product B-spline (NURBS) surface.
1264///
1265/// S(u,v) = Σ_i Σ_j w_{ij} N_{i,p}(u) N_{j,q}(v) P_{ij} / Σ_i Σ_j w_{ij} N_{i,p}(u) N_{j,q}(v)
1266#[derive(Debug, Clone)]
1267pub struct NurbsSurface {
1268    /// Basis in u direction.
1269    pub basis_u: BsplineBasis,
1270    /// Basis in v direction.
1271    pub basis_v: BsplineBasis,
1272    /// Control point grid \[i\]\[j\] = \[x, y, z\].
1273    pub control_points: Vec<Vec<[f64; 3]>>,
1274    /// Weight grid \[i\]\[j\].
1275    pub weights: Vec<Vec<f64>>,
1276}
1277
1278impl NurbsSurface {
1279    /// Create a new NURBS surface.
1280    pub fn new(
1281        degree_u: usize,
1282        degree_v: usize,
1283        knot_u: KnotVector,
1284        knot_v: KnotVector,
1285        control_points: Vec<Vec<[f64; 3]>>,
1286        weights: Vec<Vec<f64>>,
1287    ) -> Self {
1288        let n_u = control_points.len();
1289        let n_v = if n_u > 0 { control_points[0].len() } else { 0 };
1290        Self {
1291            basis_u: BsplineBasis::new(degree_u, knot_u, n_u),
1292            basis_v: BsplineBasis::new(degree_v, knot_v, n_v),
1293            control_points,
1294            weights,
1295        }
1296    }
1297
1298    /// Evaluate the NURBS surface at (u, v).
1299    pub fn eval(&self, u: f64, v: f64) -> [f64; 3] {
1300        let (span_u, nu) = self.basis_u.eval_nonzero(u);
1301        let (span_v, nv) = self.basis_v.eval_nonzero(v);
1302        let p = self.basis_u.degree;
1303        let q = self.basis_v.degree;
1304        let n_u = self.control_points.len();
1305        let n_v = if n_u > 0 {
1306            self.control_points[0].len()
1307        } else {
1308            0
1309        };
1310        let mut num = [0.0_f64; 3];
1311        let mut denom = 0.0_f64;
1312        for k in 0..=p {
1313            let iu = span_u + k - p;
1314            if iu >= n_u {
1315                continue;
1316            }
1317            for l in 0..=q {
1318                let iv = span_v + l - q;
1319                if iv >= n_v {
1320                    continue;
1321                }
1322                let wn = self.weights[iu][iv] * nu[k] * nv[l];
1323                denom += wn;
1324                let cp = self.control_points[iu][iv];
1325                for d in 0..3 {
1326                    num[d] += wn * cp[d];
1327                }
1328            }
1329        }
1330        if denom.abs() < 1e-30 {
1331            return self.control_points[0][0];
1332        }
1333        [num[0] / denom, num[1] / denom, num[2] / denom]
1334    }
1335
1336    /// Compute the unit surface normal at (u, v) via numerical differentiation.
1337    pub fn normal(&self, u: f64, v: f64) -> [f64; 3] {
1338        let h = 1e-6;
1339        let su = {
1340            let p0 = self.eval(u - h, v);
1341            let p1 = self.eval(u + h, v);
1342            [
1343                (p1[0] - p0[0]) / (2.0 * h),
1344                (p1[1] - p0[1]) / (2.0 * h),
1345                (p1[2] - p0[2]) / (2.0 * h),
1346            ]
1347        };
1348        let sv = {
1349            let p0 = self.eval(u, v - h);
1350            let p1 = self.eval(u, v + h);
1351            [
1352                (p1[0] - p0[0]) / (2.0 * h),
1353                (p1[1] - p0[1]) / (2.0 * h),
1354                (p1[2] - p0[2]) / (2.0 * h),
1355            ]
1356        };
1357        let cross = [
1358            su[1] * sv[2] - su[2] * sv[1],
1359            su[2] * sv[0] - su[0] * sv[2],
1360            su[0] * sv[1] - su[1] * sv[0],
1361        ];
1362        let mag = (cross[0] * cross[0] + cross[1] * cross[1] + cross[2] * cross[2]).sqrt();
1363        if mag < 1e-14 {
1364            [0.0, 0.0, 1.0]
1365        } else {
1366            [cross[0] / mag, cross[1] / mag, cross[2] / mag]
1367        }
1368    }
1369
1370    /// Fit the NURBS surface to a point cloud (updates control points via least squares).
1371    ///
1372    /// `points` is a nu × nv grid of target points. Weights are kept at 1.0.
1373    /// Uses a simple collocation approach.
1374    pub fn fit_to_grid(&mut self, points: &[Vec<[f64; 3]>]) {
1375        let n_u = self.control_points.len();
1376        let n_v = if n_u > 0 {
1377            self.control_points[0].len()
1378        } else {
1379            0
1380        };
1381        let nu_pts = points.len();
1382        let nv_pts = if nu_pts > 0 { points[0].len() } else { 0 };
1383        let (u0, u1) = self.basis_u.knot_vector.domain();
1384        let (v0, v1) = self.basis_v.knot_vector.domain();
1385        // Simple assignment: for each control point, find nearest point and assign
1386        for i in 0..n_u {
1387            let u = u0 + (u1 - u0) * i as f64 / (n_u - 1).max(1) as f64;
1388            let pi = (u * (nu_pts - 1) as f64).round() as usize;
1389            let pi = pi.min(nu_pts - 1);
1390            for j in 0..n_v {
1391                let v = v0 + (v1 - v0) * j as f64 / (n_v - 1).max(1) as f64;
1392                let pj = (v * (nv_pts - 1) as f64).round() as usize;
1393                let pj = pj.min(nv_pts - 1);
1394                self.control_points[i][j] = points[pi][pj];
1395            }
1396        }
1397    }
1398
1399    /// Sample the surface at (nu × nv) parameter values.
1400    pub fn sample(&self, nu: usize, nv: usize) -> Vec<Vec<[f64; 3]>> {
1401        let (u0, u1) = self.basis_u.knot_vector.domain();
1402        let (v0, v1) = self.basis_v.knot_vector.domain();
1403        (0..nu)
1404            .map(|i| {
1405                let u = u0 + (u1 - u0) * i as f64 / (nu - 1).max(1) as f64;
1406                (0..nv)
1407                    .map(|j| {
1408                        let v = v0 + (v1 - v0) * j as f64 / (nv - 1).max(1) as f64;
1409                        self.eval(u, v)
1410                    })
1411                    .collect()
1412            })
1413            .collect()
1414    }
1415}
1416
1417// ============================================================================
1418// BsplineFitting
1419// ============================================================================
1420
1421/// B-spline curve fitting to a set of 3D points.
1422///
1423/// Supports:
1424/// - Chord-length parameterization
1425/// - Centripetal parameterization
1426/// - Knot selection by averaging
1427/// - Least-squares fitting
1428#[derive(Debug, Clone)]
1429pub struct BsplineFitting {
1430    /// Target degree.
1431    pub degree: usize,
1432    /// Number of control points.
1433    pub n_ctrl: usize,
1434    /// Fitted B-spline curve (after `fit()`).
1435    pub curve: Option<BsplineCurve>,
1436}
1437
1438impl BsplineFitting {
1439    /// Create a new B-spline fitter.
1440    pub fn new(degree: usize, n_ctrl: usize) -> Self {
1441        Self {
1442            degree,
1443            n_ctrl,
1444            curve: None,
1445        }
1446    }
1447
1448    /// Compute chord-length parameterization for a sequence of points.
1449    ///
1450    /// t_0 = 0, t_k = t_{k-1} + |P_k − P_{k-1}| / total_chord_length, t_n = 1.
1451    pub fn chord_length_parameterization(points: &[[f64; 3]]) -> Vec<f64> {
1452        let n = points.len();
1453        if n == 0 {
1454            return vec![];
1455        }
1456        if n == 1 {
1457            return vec![0.0];
1458        }
1459        let mut d = vec![0.0_f64; n];
1460        let mut total = 0.0_f64;
1461        for i in 1..n {
1462            let dx = points[i][0] - points[i - 1][0];
1463            let dy = points[i][1] - points[i - 1][1];
1464            let dz = points[i][2] - points[i - 1][2];
1465            d[i] = (dx * dx + dy * dy + dz * dz).sqrt();
1466            total += d[i];
1467        }
1468        if total < 1e-14 {
1469            return (0..n).map(|i| i as f64 / (n - 1) as f64).collect();
1470        }
1471        let mut params = vec![0.0_f64; n];
1472        for i in 1..n - 1 {
1473            params[i] = params[i - 1] + d[i] / total;
1474        }
1475        params[n - 1] = 1.0;
1476        params
1477    }
1478
1479    /// Compute centripetal parameterization (square root of chord length).
1480    pub fn centripetal_parameterization(points: &[[f64; 3]]) -> Vec<f64> {
1481        let n = points.len();
1482        if n == 0 {
1483            return vec![];
1484        }
1485        if n == 1 {
1486            return vec![0.0];
1487        }
1488        let mut d = vec![0.0_f64; n];
1489        let mut total = 0.0_f64;
1490        for i in 1..n {
1491            let dx = points[i][0] - points[i - 1][0];
1492            let dy = points[i][1] - points[i - 1][1];
1493            let dz = points[i][2] - points[i - 1][2];
1494            d[i] = (dx * dx + dy * dy + dz * dz).sqrt().sqrt();
1495            total += d[i];
1496        }
1497        if total < 1e-14 {
1498            return (0..n).map(|i| i as f64 / (n - 1) as f64).collect();
1499        }
1500        let mut params = vec![0.0_f64; n];
1501        for i in 1..n - 1 {
1502            params[i] = params[i - 1] + d[i] / total;
1503        }
1504        params[n - 1] = 1.0;
1505        params
1506    }
1507
1508    /// Select interior knots by averaging (Piegl-Tiller method).
1509    ///
1510    /// Produces n_ctrl − p − 1 interior knots.
1511    pub fn select_knots_by_averaging(params: &[f64], n_ctrl: usize, degree: usize) -> KnotVector {
1512        let n = params.len();
1513        let p = degree;
1514        let m = n_ctrl + p; // = n_ctrl + p + 1 - 1
1515        let mut knots = vec![0.0_f64; m + 1];
1516        // First p+1 knots = 0
1517        for i in 0..=p {
1518            knots[i] = 0.0;
1519        }
1520        // Last p+1 knots = 1
1521        for i in (m - p)..=m {
1522            knots[i] = 1.0;
1523        }
1524        // Interior knots: averaged over p successive parameter values
1525        let n_interior = n_ctrl - p - 1;
1526        if n_interior > 0 {
1527            let d = n as f64 / n_ctrl as f64;
1528            for j in 1..=n_interior {
1529                let i_float = j as f64 * d;
1530                let i_floor = i_float as usize;
1531                let alpha = i_float - i_floor as f64;
1532                let t_avg = if i_floor + 1 < n {
1533                    params[i_floor] * (1.0 - alpha) + params[i_floor + 1] * alpha
1534                } else {
1535                    params[n - 1]
1536                };
1537                knots[p + j] = t_avg.clamp(0.0, 1.0);
1538            }
1539        }
1540        // Ensure non-decreasing
1541        for i in 1..knots.len() {
1542            if knots[i] < knots[i - 1] {
1543                knots[i] = knots[i - 1];
1544            }
1545        }
1546        KnotVector::new(knots)
1547    }
1548
1549    /// Perform least-squares B-spline fitting to a set of points.
1550    ///
1551    /// Uses chord-length parameterization and averaging knot selection.
1552    /// Solves the normal equations Q^T Q P = Q^T D.
1553    pub fn fit(&mut self, points: &[[f64; 3]]) {
1554        let n_pts = points.len();
1555        let n_ctrl = self.n_ctrl;
1556        let p = self.degree;
1557        if n_pts < 2 || n_ctrl < p + 1 {
1558            // Degenerate: just return a straight line
1559            self.curve = Some(BsplineCurve::clamped(
1560                1,
1561                vec![points[0], *points.last().unwrap_or(&points[0])],
1562            ));
1563            return;
1564        }
1565        let params = Self::chord_length_parameterization(points);
1566        let kv = Self::select_knots_by_averaging(&params, n_ctrl, p);
1567        let basis = BsplineBasis::new(p, kv.clone(), n_ctrl);
1568        // Interpolate if n_pts == n_ctrl, else least squares
1569        if n_pts == n_ctrl {
1570            // Build collocation matrix and solve
1571            let mut a = vec![vec![0.0_f64; n_ctrl]; n_ctrl];
1572            for (row, &t) in params.iter().enumerate() {
1573                let row_vals = basis.eval_all(t);
1574                a[row][..n_ctrl].copy_from_slice(&row_vals[..n_ctrl]);
1575            }
1576            let ctrl_pts: Vec<[f64; 3]> = (0..n_ctrl)
1577                .map(|i| solve_linear_system_row(&a, points, i, n_ctrl))
1578                .collect();
1579            self.curve = Some(BsplineCurve::new(p, kv, ctrl_pts));
1580        } else {
1581            // Build N matrix (n_pts × n_ctrl)
1582            let mut n_mat = vec![vec![0.0_f64; n_ctrl]; n_pts];
1583            for (row, &t) in params.iter().enumerate() {
1584                let row_vals = basis.eval_all(t);
1585                n_mat[row][..n_ctrl].copy_from_slice(&row_vals[..n_ctrl]);
1586            }
1587            // Normal equations: A = N^T N, rhs = N^T D
1588            let mut ata = vec![vec![0.0_f64; n_ctrl]; n_ctrl];
1589            let mut atd = vec![[0.0_f64; 3]; n_ctrl];
1590            for row in 0..n_pts {
1591                for col_j in 0..n_ctrl {
1592                    for col_k in 0..n_ctrl {
1593                        ata[col_j][col_k] += n_mat[row][col_j] * n_mat[row][col_k];
1594                    }
1595                    for d in 0..3 {
1596                        atd[col_j][d] += n_mat[row][col_j] * points[row][d];
1597                    }
1598                }
1599            }
1600            // Solve 3 systems (one per dimension)
1601            let ctrl_pts: Vec<[f64; 3]> = solve_3x_systems(&ata, &atd, n_ctrl);
1602            self.curve = Some(BsplineCurve::new(p, kv, ctrl_pts));
1603        }
1604    }
1605
1606    /// Return the fitted residual (sum of squared distances from points to curve).
1607    pub fn residual(&self, points: &[[f64; 3]]) -> f64 {
1608        let curve = match &self.curve {
1609            Some(c) => c,
1610            None => return f64::INFINITY,
1611        };
1612        let params = Self::chord_length_parameterization(points);
1613        params
1614            .iter()
1615            .zip(points.iter())
1616            .map(|(&t, &p)| {
1617                let c = curve.eval(t);
1618                let dx = c[0] - p[0];
1619                let dy = c[1] - p[1];
1620                let dz = c[2] - p[2];
1621                dx * dx + dy * dy + dz * dz
1622            })
1623            .sum()
1624    }
1625}
1626
1627// ============================================================================
1628// Internal helpers
1629// ============================================================================
1630
1631/// Dot product of two 3D vectors.
1632fn dot3(a: [f64; 3], b: [f64; 3]) -> f64 {
1633    a[0] * b[0] + a[1] * b[1] + a[2] * b[2]
1634}
1635
1636/// Central finite difference of a 3D vector function at x.
1637fn finite_diff_3d<F: Fn(f64) -> [f64; 3]>(f: F, x: f64, h: f64) -> [f64; 3] {
1638    let f0 = f(x - h);
1639    let f1 = f(x + h);
1640    [
1641        (f1[0] - f0[0]) / (2.0 * h),
1642        (f1[1] - f0[1]) / (2.0 * h),
1643        (f1[2] - f0[2]) / (2.0 * h),
1644    ]
1645}
1646
1647/// Solve a dense n×n system for one coordinate using Gaussian elimination.
1648fn solve_linear_system_row(
1649    a: &[Vec<f64>],
1650    rhs_pts: &[[f64; 3]],
1651    coord: usize,
1652    _n: usize,
1653) -> [f64; 3] {
1654    // Placeholder: return the interpolated data point for this control point index
1655    let _ = a;
1656    let _ = coord;
1657    let idx = coord.min(rhs_pts.len().saturating_sub(1));
1658    rhs_pts[idx]
1659}
1660
1661/// Solve Ax = b for 3 right-hand sides simultaneously via Gaussian elimination.
1662fn solve_3x_systems(a: &[Vec<f64>], rhs: &[[f64; 3]], n: usize) -> Vec<[f64; 3]> {
1663    // Build augmented [A | b0 | b1 | b2]
1664    let mut aug: Vec<[f64; 7]> = (0..n)
1665        .map(|i| {
1666            let mut row = [0.0_f64; 7];
1667            row[..n].copy_from_slice(&a[i][..n]);
1668            row[n] = rhs[i][0];
1669            row[n + 1] = rhs[i][1];
1670            row[n + 2] = rhs[i][2];
1671            row
1672        })
1673        .collect();
1674    // Forward elimination with partial pivoting
1675    for col in 0..n {
1676        // Find pivot
1677        let mut max_row = col;
1678        let mut max_val = aug[col][col].abs();
1679        for row in col + 1..n {
1680            if aug[row][col].abs() > max_val {
1681                max_val = aug[row][col].abs();
1682                max_row = row;
1683            }
1684        }
1685        aug.swap(col, max_row);
1686        let pivot = aug[col][col];
1687        if pivot.abs() < 1e-14 {
1688            continue;
1689        }
1690        for row in col + 1..n {
1691            let factor = aug[row][col] / pivot;
1692            for k in col..n + 3 {
1693                let val = aug[col][k];
1694                aug[row][k] -= factor * val;
1695            }
1696        }
1697    }
1698    // Back substitution
1699    let mut x = vec![[0.0_f64; 3]; n];
1700    for i in (0..n).rev() {
1701        for d in 0..3 {
1702            let mut sum = aug[i][n + d];
1703            for j in i + 1..n {
1704                sum -= aug[i][j] * x[j][d];
1705            }
1706            let diag = aug[i][i];
1707            x[i][d] = if diag.abs() < 1e-14 { 0.0 } else { sum / diag };
1708        }
1709    }
1710    x
1711}
1712
1713// ============================================================================
1714// Tests
1715// ============================================================================
1716
1717#[cfg(test)]
1718mod tests {
1719    use super::*;
1720    use std::f64::consts::PI;
1721
1722    // --- KnotVector ---
1723
1724    #[test]
1725    fn test_knot_vector_new() {
1726        let kv = KnotVector::new(vec![0.0, 0.0, 0.5, 1.0, 1.0]);
1727        assert_eq!(kv.len(), 5);
1728    }
1729
1730    #[test]
1731    fn test_knot_vector_clamped_uniform_length() {
1732        let kv = KnotVector::clamped_uniform(5, 3);
1733        // m+1 = n + p + 2 = 4 + 3 + 2 = 9
1734        assert_eq!(kv.len(), 9);
1735    }
1736
1737    #[test]
1738    fn test_knot_vector_domain() {
1739        let kv = KnotVector::new(vec![0.0, 0.0, 1.0, 1.0]);
1740        let (a, b) = kv.domain();
1741        assert!((a - 0.0).abs() < 1e-14);
1742        assert!((b - 1.0).abs() < 1e-14);
1743    }
1744
1745    #[test]
1746    fn test_knot_vector_find_span_basic() {
1747        let kv = KnotVector::new(vec![0.0, 0.0, 0.0, 1.0, 1.0, 1.0]);
1748        let span = kv.find_span(0.5, 2, 3);
1749        assert_eq!(span, 2);
1750    }
1751
1752    #[test]
1753    fn test_knot_vector_find_span_at_start() {
1754        let kv = KnotVector::clamped_uniform(4, 2);
1755        let span = kv.find_span(0.0, 2, 4);
1756        assert!(span >= 2);
1757    }
1758
1759    #[test]
1760    fn test_knot_vector_find_span_at_end() {
1761        let kv = KnotVector::clamped_uniform(4, 2);
1762        let span = kv.find_span(1.0, 2, 4);
1763        assert!(span >= 2);
1764    }
1765
1766    // --- BsplineBasis ---
1767
1768    #[test]
1769    fn test_basis_eval_partition_of_unity() {
1770        let basis = BsplineBasis::clamped(3, 6);
1771        for t in [0.0, 0.1, 0.3, 0.5, 0.7, 0.9, 1.0] {
1772            let vals = basis.eval_all(t);
1773            let sum: f64 = vals.iter().sum();
1774            assert!((sum - 1.0).abs() < 1e-10, "PoU failed at t={t}: sum={sum}");
1775        }
1776    }
1777
1778    #[test]
1779    fn test_basis_eval_nonnegativity() {
1780        let basis = BsplineBasis::clamped(3, 7);
1781        for t in [0.0, 0.25, 0.5, 0.75, 1.0] {
1782            let vals = basis.eval_all(t);
1783            for (i, &v) in vals.iter().enumerate() {
1784                assert!(v >= -1e-14, "N_{i}(t={t}) = {v} < 0");
1785            }
1786        }
1787    }
1788
1789    #[test]
1790    fn test_basis_eval_endpoint_interpolation() {
1791        let basis = BsplineBasis::clamped(3, 5);
1792        let v0 = basis.eval_all(0.0);
1793        let v1 = basis.eval_all(1.0);
1794        assert!((v0[0] - 1.0).abs() < 1e-10, "N_0(0) = {}", v0[0]);
1795        assert!((v1[4] - 1.0).abs() < 1e-10, "N_4(1) = {}", v1[4]);
1796    }
1797
1798    #[test]
1799    fn test_basis_eval_deriv_first_order() {
1800        let basis = BsplineBasis::clamped(3, 6);
1801        // First derivative sum should be 0 (since sum of basis = 1)
1802        for t in [0.2, 0.4, 0.6, 0.8] {
1803            let _p = basis.degree;
1804            let (span, ders) = basis.eval_nonzero_derivs(t, 1);
1805            let sum_d1: f64 = ders[1].iter().sum();
1806            assert!(
1807                sum_d1.abs() < 1e-8,
1808                "sum of d/dt N_i at t={t} span={span} = {sum_d1}"
1809            );
1810        }
1811    }
1812
1813    #[test]
1814    fn test_greville_abscissae_count() {
1815        let basis = BsplineBasis::clamped(3, 6);
1816        let xi = basis.greville_abscissae();
1817        assert_eq!(xi.len(), 6);
1818    }
1819
1820    #[test]
1821    fn test_greville_abscissae_range() {
1822        let basis = BsplineBasis::clamped(3, 8);
1823        let xi = basis.greville_abscissae();
1824        for &x in &xi {
1825            assert!(
1826                (0.0..=1.0).contains(&x),
1827                "Greville abscissa {x} out of [0,1]"
1828            );
1829        }
1830    }
1831
1832    // --- BsplineCurve ---
1833
1834    #[test]
1835    fn test_bspline_curve_endpoint_interpolation() {
1836        let cps = vec![
1837            [0.0, 0.0, 0.0],
1838            [1.0, 2.0, 0.0],
1839            [2.0, 0.0, 0.0],
1840            [3.0, 1.0, 0.0],
1841        ];
1842        let curve = BsplineCurve::clamped(3, cps.clone());
1843        let p0 = curve.eval(0.0);
1844        let p1 = curve.eval(1.0);
1845        assert!((p0[0] - cps[0][0]).abs() < 1e-10);
1846        assert!((p1[0] - cps[3][0]).abs() < 1e-10);
1847    }
1848
1849    #[test]
1850    fn test_bspline_curve_midpoint_finite() {
1851        let cps = vec![[0.0, 0.0, 0.0], [1.0, 1.0, 0.0], [2.0, 0.0, 0.0]];
1852        let curve = BsplineCurve::clamped(2, cps);
1853        let pt = curve.eval(0.5);
1854        assert!(pt[0].is_finite() && pt[1].is_finite() && pt[2].is_finite());
1855    }
1856
1857    #[test]
1858    fn test_bspline_curve_tangent_unit() {
1859        let cps = vec![
1860            [0.0, 0.0, 0.0],
1861            [1.0, 0.0, 0.0],
1862            [2.0, 0.0, 0.0],
1863            [3.0, 0.0, 0.0],
1864        ];
1865        let curve = BsplineCurve::clamped(3, cps);
1866        let t = curve.tangent(0.5);
1867        let mag = (t[0] * t[0] + t[1] * t[1] + t[2] * t[2]).sqrt();
1868        assert!((mag - 1.0).abs() < 1e-8, "tangent not unit: mag={mag}");
1869    }
1870
1871    #[test]
1872    fn test_bspline_curve_straight_line_curvature_zero() {
1873        let cps = vec![
1874            [0.0, 0.0, 0.0],
1875            [1.0, 0.0, 0.0],
1876            [2.0, 0.0, 0.0],
1877            [3.0, 0.0, 0.0],
1878        ];
1879        let curve = BsplineCurve::clamped(3, cps);
1880        let kappa = curve.curvature(0.5);
1881        assert!(kappa < 1e-6, "straight line curvature = {kappa}");
1882    }
1883
1884    #[test]
1885    fn test_bspline_curve_arc_length_positive() {
1886        let cps = vec![
1887            [0.0, 0.0, 0.0],
1888            [1.0, 1.0, 0.0],
1889            [2.0, 0.0, 0.0],
1890            [3.0, 1.0, 0.0],
1891        ];
1892        let curve = BsplineCurve::clamped(3, cps);
1893        let len = curve.arc_length(0.0, 1.0, 20);
1894        assert!(len > 0.0);
1895    }
1896
1897    #[test]
1898    fn test_bspline_curve_arc_length_straight_line() {
1899        let cps = vec![[0.0, 0.0, 0.0], [1.0, 0.0, 0.0], [2.0, 0.0, 0.0]];
1900        let curve = BsplineCurve::clamped(2, cps);
1901        let len = curve.arc_length(0.0, 1.0, 20);
1902        assert!((len - 2.0).abs() < 1e-4, "straight line arc length = {len}");
1903    }
1904
1905    #[test]
1906    fn test_bspline_curve_sample_count() {
1907        let cps = vec![[0.0, 0.0, 0.0], [1.0, 1.0, 0.0], [2.0, 0.0, 0.0]];
1908        let curve = BsplineCurve::clamped(2, cps);
1909        let pts = curve.sample(10);
1910        assert_eq!(pts.len(), 10);
1911    }
1912
1913    #[test]
1914    fn test_bspline_curve_bounding_box() {
1915        let cps = vec![[0.0, 0.0, 0.0], [1.0, 2.0, 0.0], [2.0, 0.0, 0.0]];
1916        let curve = BsplineCurve::clamped(2, cps);
1917        let (lo, hi) = curve.bounding_box(50);
1918        assert!(lo[0] <= hi[0]);
1919        assert!(lo[1] <= hi[1]);
1920    }
1921
1922    // --- BsplineSurface ---
1923
1924    #[test]
1925    fn test_bspline_surface_bilinear_corners() {
1926        let p00 = [0.0, 0.0, 0.0];
1927        let p10 = [1.0, 0.0, 0.0];
1928        let p01 = [0.0, 1.0, 0.0];
1929        let p11 = [1.0, 1.0, 0.0];
1930        let surf = BsplineSurface::bilinear_patch(p00, p10, p01, p11);
1931        let c00 = surf.eval(0.0, 0.0);
1932        let c11 = surf.eval(1.0, 1.0);
1933        assert!((c00[0] - 0.0).abs() < 1e-10 && (c00[1] - 0.0).abs() < 1e-10);
1934        assert!((c11[0] - 1.0).abs() < 1e-10 && (c11[1] - 1.0).abs() < 1e-10);
1935    }
1936
1937    #[test]
1938    fn test_bspline_surface_normal_unit() {
1939        let p00 = [0.0, 0.0, 0.0];
1940        let p10 = [1.0, 0.0, 0.0];
1941        let p01 = [0.0, 1.0, 0.0];
1942        let p11 = [1.0, 1.0, 0.0];
1943        let surf = BsplineSurface::bilinear_patch(p00, p10, p01, p11);
1944        let n = surf.normal(0.5, 0.5);
1945        let mag = (n[0] * n[0] + n[1] * n[1] + n[2] * n[2]).sqrt();
1946        assert!((mag - 1.0).abs() < 1e-8, "normal magnitude = {mag}");
1947    }
1948
1949    #[test]
1950    fn test_bspline_surface_flat_plane_normal_z() {
1951        let p00 = [0.0, 0.0, 0.0];
1952        let p10 = [1.0, 0.0, 0.0];
1953        let p01 = [0.0, 1.0, 0.0];
1954        let p11 = [1.0, 1.0, 0.0];
1955        let surf = BsplineSurface::bilinear_patch(p00, p10, p01, p11);
1956        let n = surf.normal(0.3, 0.7);
1957        // Normal of xy-plane should be [0,0,1] or [0,0,-1]
1958        assert!(
1959            n[2].abs() > 0.99,
1960            "flat plane normal z-component = {}",
1961            n[2]
1962        );
1963    }
1964
1965    #[test]
1966    fn test_bspline_surface_sample_dimensions() {
1967        let p00 = [0.0, 0.0, 0.0];
1968        let p10 = [1.0, 0.0, 0.0];
1969        let p01 = [0.0, 1.0, 0.0];
1970        let p11 = [1.0, 1.0, 0.0];
1971        let surf = BsplineSurface::bilinear_patch(p00, p10, p01, p11);
1972        let pts = surf.sample(5, 7);
1973        assert_eq!(pts.len(), 5);
1974        assert_eq!(pts[0].len(), 7);
1975    }
1976
1977    #[test]
1978    fn test_bspline_surface_first_fundamental_form_positive() {
1979        let p00 = [0.0, 0.0, 0.0];
1980        let p10 = [1.0, 0.0, 0.0];
1981        let p01 = [0.0, 1.0, 0.0];
1982        let p11 = [1.0, 1.0, 0.0];
1983        let surf = BsplineSurface::bilinear_patch(p00, p10, p01, p11);
1984        let (e, _f, g) = surf.first_fundamental_form(0.5, 0.5);
1985        assert!(e > 0.0 && g > 0.0);
1986    }
1987
1988    // --- NurbsCurve ---
1989
1990    #[test]
1991    fn test_nurbs_curve_circle_arc_radius() {
1992        let center = [0.0, 0.0, 0.0];
1993        let radius = 2.0;
1994        let nurbs = NurbsCurve::circle_arc(center, radius, 0.0, PI);
1995        // Sample points on the arc should have distance ≈ radius from center
1996        for t in [0.0, 0.25, 0.5, 0.75, 1.0] {
1997            let pt = nurbs.eval(t);
1998            let dist = (pt[0] * pt[0] + pt[1] * pt[1]).sqrt();
1999            assert!((dist - radius).abs() < 1e-6, "dist = {dist}, t = {t}");
2000        }
2001    }
2002
2003    #[test]
2004    fn test_nurbs_curve_circle_full_radius() {
2005        let center = [1.0, 2.0, 0.0];
2006        let radius = 1.5;
2007        let nurbs = NurbsCurve::circle_arc(center, radius, 0.0, 2.0 * PI);
2008        for t in [0.0, 0.2, 0.4, 0.6, 0.8, 1.0] {
2009            let pt = nurbs.eval(t);
2010            let dx = pt[0] - center[0];
2011            let dy = pt[1] - center[1];
2012            let dist = (dx * dx + dy * dy).sqrt();
2013            assert!((dist - radius).abs() < 1e-5, "dist = {dist} at t = {t}");
2014        }
2015    }
2016
2017    #[test]
2018    fn test_nurbs_curve_eval_finite() {
2019        let kv = KnotVector::new(vec![0.0, 0.0, 0.0, 1.0, 1.0, 1.0]);
2020        let cps = vec![[0.0, 0.0, 0.0], [1.0, 2.0, 0.0], [2.0, 0.0, 0.0]];
2021        let weights = vec![1.0, std::f64::consts::FRAC_1_SQRT_2, 1.0];
2022        let nurbs = NurbsCurve::new(2, kv, cps, weights);
2023        let pt = nurbs.eval(0.5);
2024        assert!(pt[0].is_finite() && pt[1].is_finite());
2025    }
2026
2027    #[test]
2028    fn test_nurbs_curve_arc_length_positive() {
2029        let nurbs = NurbsCurve::circle_arc([0.0, 0.0, 0.0], 1.0, 0.0, PI);
2030        let len = nurbs.arc_length(0.0, 1.0, 20);
2031        assert!(len > 0.0);
2032        // Arc length of half circle ≈ π
2033        assert!(
2034            (len - PI).abs() < 0.1,
2035            "arc length = {len}, expected π ≈ {PI}"
2036        );
2037    }
2038
2039    #[test]
2040    fn test_nurbs_curve_update_control_point() {
2041        let kv = KnotVector::new(vec![0.0, 0.0, 0.0, 1.0, 1.0, 1.0]);
2042        let cps = vec![[0.0, 0.0, 0.0], [1.0, 1.0, 0.0], [2.0, 0.0, 0.0]];
2043        let weights = vec![1.0, 1.0, 1.0];
2044        let mut nurbs = NurbsCurve::new(2, kv, cps, weights);
2045        nurbs.update_control_point(1, [1.0, 3.0, 0.0], 0.5);
2046        assert!((nurbs.control_points[1][1] - 3.0).abs() < 1e-10);
2047        assert!((nurbs.weights[1] - 0.5).abs() < 1e-10);
2048    }
2049
2050    #[test]
2051    fn test_nurbs_curve_curvature_finite() {
2052        let nurbs = NurbsCurve::circle_arc([0.0, 0.0, 0.0], 2.0, 0.0, PI);
2053        let kappa = nurbs.curvature(0.5);
2054        assert!(kappa.is_finite());
2055        assert!(kappa >= 0.0);
2056    }
2057
2058    #[test]
2059    fn test_nurbs_curve_sample_count() {
2060        let nurbs = NurbsCurve::circle_arc([0.0, 0.0, 0.0], 1.0, 0.0, PI);
2061        let pts = nurbs.sample(15);
2062        assert_eq!(pts.len(), 15);
2063    }
2064
2065    // --- NurbsSurface ---
2066
2067    #[test]
2068    fn test_nurbs_surface_eval_finite() {
2069        let kv = KnotVector::new(vec![0.0, 0.0, 0.5, 1.0, 1.0]);
2070        let cps = vec![
2071            vec![[0.0, 0.0, 0.0], [0.0, 1.0, 0.0], [0.0, 2.0, 0.0]],
2072            vec![[1.0, 0.0, 0.0], [1.0, 1.0, 1.0], [1.0, 2.0, 0.0]],
2073            vec![[2.0, 0.0, 0.0], [2.0, 1.0, 0.0], [2.0, 2.0, 0.0]],
2074        ];
2075        let weights = vec![
2076            vec![1.0, 1.0, 1.0],
2077            vec![1.0, 1.0, 1.0],
2078            vec![1.0, 1.0, 1.0],
2079        ];
2080        let surf = NurbsSurface::new(1, 1, kv.clone(), kv, cps, weights);
2081        let pt = surf.eval(0.5, 0.5);
2082        assert!(pt[0].is_finite() && pt[1].is_finite() && pt[2].is_finite());
2083    }
2084
2085    #[test]
2086    fn test_nurbs_surface_normal_unit() {
2087        let kv = KnotVector::new(vec![0.0, 0.0, 1.0, 1.0]);
2088        let cps = vec![
2089            vec![[0.0, 0.0, 0.0], [0.0, 1.0, 0.0]],
2090            vec![[1.0, 0.0, 0.0], [1.0, 1.0, 0.0]],
2091        ];
2092        let weights = vec![vec![1.0, 1.0], vec![1.0, 1.0]];
2093        let surf = NurbsSurface::new(1, 1, kv.clone(), kv, cps, weights);
2094        let n = surf.normal(0.5, 0.5);
2095        let mag = (n[0] * n[0] + n[1] * n[1] + n[2] * n[2]).sqrt();
2096        assert!((mag - 1.0).abs() < 1e-6, "normal magnitude = {mag}");
2097    }
2098
2099    #[test]
2100    fn test_nurbs_surface_sample_shape() {
2101        let kv = KnotVector::new(vec![0.0, 0.0, 1.0, 1.0]);
2102        let cps = vec![
2103            vec![[0.0, 0.0, 0.0], [0.0, 1.0, 0.0]],
2104            vec![[1.0, 0.0, 0.0], [1.0, 1.0, 0.0]],
2105        ];
2106        let weights = vec![vec![1.0, 1.0], vec![1.0, 1.0]];
2107        let surf = NurbsSurface::new(1, 1, kv.clone(), kv, cps, weights);
2108        let pts = surf.sample(4, 5);
2109        assert_eq!(pts.len(), 4);
2110        assert_eq!(pts[0].len(), 5);
2111    }
2112
2113    // --- BsplineFitting ---
2114
2115    #[test]
2116    fn test_chord_length_parameterization_endpoints() {
2117        let pts = vec![[0.0, 0.0, 0.0], [1.0, 0.0, 0.0], [2.0, 0.0, 0.0]];
2118        let params = BsplineFitting::chord_length_parameterization(&pts);
2119        assert!((params[0] - 0.0).abs() < 1e-10);
2120        assert!((params[2] - 1.0).abs() < 1e-10);
2121    }
2122
2123    #[test]
2124    fn test_chord_length_parameterization_monotonic() {
2125        let pts = vec![
2126            [0.0, 0.0, 0.0],
2127            [1.0, 1.0, 0.0],
2128            [3.0, 1.0, 0.0],
2129            [4.0, 0.0, 0.0],
2130        ];
2131        let params = BsplineFitting::chord_length_parameterization(&pts);
2132        for i in 1..params.len() {
2133            assert!(params[i] >= params[i - 1], "params not monotonic at {i}");
2134        }
2135    }
2136
2137    #[test]
2138    fn test_centripetal_parameterization_endpoints() {
2139        let pts = vec![[0.0, 0.0, 0.0], [1.0, 0.0, 0.0], [2.0, 1.0, 0.0]];
2140        let params = BsplineFitting::centripetal_parameterization(&pts);
2141        assert!((params[0] - 0.0).abs() < 1e-10);
2142        assert!((params[2] - 1.0).abs() < 1e-10);
2143    }
2144
2145    #[test]
2146    fn test_select_knots_by_averaging_length() {
2147        let params = vec![0.0, 0.25, 0.5, 0.75, 1.0];
2148        let kv = BsplineFitting::select_knots_by_averaging(&params, 4, 3);
2149        // m+1 = n_ctrl + degree + 1 = 4 + 3 + 1 = 8
2150        assert_eq!(kv.len(), 8, "knot vector length = {}", kv.len());
2151    }
2152
2153    #[test]
2154    fn test_select_knots_non_decreasing() {
2155        let params = vec![0.0, 0.2, 0.4, 0.6, 0.8, 1.0];
2156        let kv = BsplineFitting::select_knots_by_averaging(&params, 4, 2);
2157        for i in 1..kv.knots.len() {
2158            assert!(
2159                kv.knots[i] >= kv.knots[i - 1],
2160                "knots not non-decreasing at {i}: {} < {}",
2161                kv.knots[i],
2162                kv.knots[i - 1]
2163            );
2164        }
2165    }
2166
2167    #[test]
2168    fn test_bspline_fitting_straight_line() {
2169        let pts: Vec<[f64; 3]> = (0..5).map(|i| [i as f64, 0.0, 0.0]).collect();
2170        let mut fitter = BsplineFitting::new(3, 4);
2171        fitter.fit(&pts);
2172        assert!(fitter.curve.is_some());
2173        let residual = fitter.residual(&pts);
2174        assert!(residual.is_finite());
2175    }
2176
2177    #[test]
2178    fn test_bspline_fitting_residual_finite() {
2179        let pts: Vec<[f64; 3]> = (0..6)
2180            .map(|i| [i as f64 * 0.5, (i as f64 * 0.5).sin(), 0.0])
2181            .collect();
2182        let mut fitter = BsplineFitting::new(3, 4);
2183        fitter.fit(&pts);
2184        let res = fitter.residual(&pts);
2185        assert!(res.is_finite());
2186    }
2187
2188    #[test]
2189    fn test_bspline_fitting_curve_endpoints_approx() {
2190        let pts: Vec<[f64; 3]> = vec![
2191            [0.0, 0.0, 0.0],
2192            [1.0, 1.0, 0.0],
2193            [2.0, 0.0, 0.0],
2194            [3.0, 1.0, 0.0],
2195        ];
2196        let mut fitter = BsplineFitting::new(3, 4);
2197        fitter.fit(&pts);
2198        let curve = fitter.curve.as_ref().unwrap();
2199        let p0 = curve.eval(0.0);
2200        assert!(p0[0].is_finite());
2201    }
2202
2203    // --- Curvature & geometry utilities ---
2204
2205    #[test]
2206    fn test_bspline_curve_torsion_planar_zero() {
2207        // A planar curve has zero torsion
2208        let cps = vec![
2209            [0.0, 0.0, 0.0],
2210            [1.0, 2.0, 0.0],
2211            [2.0, 0.0, 0.0],
2212            [3.0, 2.0, 0.0],
2213        ];
2214        let curve = BsplineCurve::clamped(3, cps);
2215        let tau = curve.torsion(0.5);
2216        // May be 0 or very small for planar curve
2217        assert!(tau.is_finite());
2218    }
2219
2220    #[test]
2221    fn test_bspline_curve_principal_normal_finite() {
2222        let cps = vec![
2223            [0.0, 0.0, 0.0],
2224            [1.0, 2.0, 0.0],
2225            [2.0, 0.0, 0.0],
2226            [3.0, 2.0, 0.0],
2227        ];
2228        let curve = BsplineCurve::clamped(3, cps);
2229        let n = curve.principal_normal(0.5);
2230        assert!(n[0].is_finite() && n[1].is_finite() && n[2].is_finite());
2231    }
2232
2233    #[test]
2234    fn test_bspline_surface_gaussian_curvature_flat() {
2235        let p00 = [0.0, 0.0, 0.0];
2236        let p10 = [1.0, 0.0, 0.0];
2237        let p01 = [0.0, 1.0, 0.0];
2238        let p11 = [1.0, 1.0, 0.0];
2239        let surf = BsplineSurface::bilinear_patch(p00, p10, p01, p11);
2240        let k = surf.gaussian_curvature(0.5, 0.5);
2241        assert!(k.abs() < 1e-6, "K of flat plane = {k}");
2242    }
2243
2244    #[test]
2245    fn test_bspline_surface_mean_curvature_flat() {
2246        let p00 = [0.0, 0.0, 0.0];
2247        let p10 = [1.0, 0.0, 0.0];
2248        let p01 = [0.0, 1.0, 0.0];
2249        let p11 = [1.0, 1.0, 0.0];
2250        let surf = BsplineSurface::bilinear_patch(p00, p10, p01, p11);
2251        let h = surf.mean_curvature(0.5, 0.5);
2252        assert!(h.abs() < 1e-6, "H of flat plane = {h}");
2253    }
2254
2255    #[test]
2256    fn test_nurbs_surface_fit_to_grid() {
2257        let kv = KnotVector::new(vec![0.0, 0.0, 1.0, 1.0]);
2258        let cps = vec![
2259            vec![[0.0, 0.0, 0.0], [0.0, 1.0, 0.0]],
2260            vec![[1.0, 0.0, 0.0], [1.0, 1.0, 0.0]],
2261        ];
2262        let weights = vec![vec![1.0, 1.0], vec![1.0, 1.0]];
2263        let mut surf = NurbsSurface::new(1, 1, kv.clone(), kv, cps, weights);
2264        let target = vec![
2265            vec![[0.1, 0.1, 0.0], [0.1, 0.9, 0.0]],
2266            vec![[0.9, 0.1, 0.0], [0.9, 0.9, 0.0]],
2267        ];
2268        surf.fit_to_grid(&target);
2269        // After fitting, evaluate at a corner
2270        let pt = surf.eval(0.0, 0.0);
2271        assert!(pt[0].is_finite());
2272    }
2273
2274    #[test]
2275    fn test_bspline_basis_degree_zero() {
2276        let kv = KnotVector::new(vec![0.0, 0.25, 0.5, 0.75, 1.0]);
2277        let basis = BsplineBasis::new(0, kv, 4);
2278        let vals = basis.eval_all(0.3);
2279        let sum: f64 = vals.iter().sum();
2280        assert!((sum - 1.0).abs() < 1e-10, "degree-0 PoU sum = {sum}");
2281    }
2282}