Skip to main content

threecrate_algorithms/
gicp.rs

1//! Generalized ICP (GICP) registration algorithm.
2//!
3//! GICP (Segal et al., 2009) extends standard ICP by modelling each point as a local
4//! Gaussian distribution whose covariance is estimated from its k nearest neighbours.
5//! Each correspondence is weighted by the combined source + target covariance, giving
6//! much stronger robustness to noise than point-to-point or point-to-plane ICP.
7//!
8//! Used in: LOAM, LIO-SAM, and most modern LiDAR odometry pipelines.
9
10use rayon::prelude::*;
11use nalgebra::{Matrix3, Matrix6, Translation3, UnitQuaternion, Vector6};
12use threecrate_core::{Error, Isometry3, NearestNeighborSearch, Point3f, PointCloud, Result, Vector3f};
13
14use crate::nearest_neighbor::KdTree;
15use crate::registration::ICPResult;
16
17// ---------------------------------------------------------------------------
18// Config
19// ---------------------------------------------------------------------------
20
21/// Parameters for [`gicp`].
22#[derive(Debug, Clone)]
23pub struct GicpConfig {
24    /// Maximum number of ICP iterations.
25    pub max_iterations: usize,
26    /// Maximum Euclidean distance for accepting a source–target correspondence.
27    pub max_correspondence_distance: f32,
28    /// Convergence threshold: stop when |ΔMSE| < this value.
29    pub convergence_threshold: f32,
30    /// Number of nearest neighbours used to estimate per-point covariance matrices.
31    pub k_correspondences: usize,
32}
33
34impl Default for GicpConfig {
35    fn default() -> Self {
36        Self {
37            max_iterations: 50,
38            max_correspondence_distance: 1.0,
39            convergence_threshold: 1e-6,
40            k_correspondences: 20,
41        }
42    }
43}
44
45// ---------------------------------------------------------------------------
46// Internal helpers
47// ---------------------------------------------------------------------------
48
49/// Skew-symmetric ("cross-product") matrix of a 3-D vector.
50#[inline]
51fn skew_sym(v: &nalgebra::Vector3<f32>) -> Matrix3<f32> {
52    Matrix3::new(
53         0.0, -v.z,  v.y,
54         v.z,  0.0, -v.x,
55        -v.y,  v.x,  0.0,
56    )
57}
58
59/// Estimate per-point 3×3 covariance matrices from k nearest neighbours.
60///
61/// Returns a `Vec` of length `points.len()`, aligned index-for-index.
62fn compute_covariances(points: &[Point3f], k: usize) -> Result<Vec<Matrix3<f32>>> {
63    let k = k.max(4); // need ≥ 4 neighbours for a non-degenerate 3-D covariance
64    let tree = KdTree::new(points)?;
65
66    let covs: Vec<Matrix3<f32>> = points
67        .par_iter()
68        .map(|p| {
69            let neighbours = tree.find_k_nearest(p, k);
70            let n = neighbours.len();
71            if n < 3 {
72                // Sparse neighbourhood: isotropic fallback.
73                return Matrix3::identity() * 1e-3_f32;
74            }
75
76            let nf = n as f32;
77            let mean = neighbours
78                .iter()
79                .map(|(idx, _)| points[*idx].coords)
80                .fold(nalgebra::Vector3::zeros(), |acc, v| acc + v)
81                / nf;
82
83            let mut cov = Matrix3::zeros();
84            for (idx, _) in &neighbours {
85                let d = points[*idx].coords - mean;
86                cov += d * d.transpose();
87            }
88            let mut cov = cov / (nf - 1.0).max(1.0);
89            // Regularise: add a small isotropic term so that locally planar or
90            // collinear neighbourhoods still have an invertible covariance matrix.
91            cov += Matrix3::identity() * 1e-4_f32;
92            cov
93        })
94        .collect();
95
96    Ok(covs)
97}
98
99// ---------------------------------------------------------------------------
100// Public API
101// ---------------------------------------------------------------------------
102
103/// Align `source` to `target` using Generalized ICP.
104///
105/// # Algorithm
106///
107/// 1. Estimate per-point covariances for both clouds using k-NN.
108/// 2. Iteratively:
109///    a. Transform source with the current estimate.
110///    b. Find nearest-neighbour correspondences within `max_correspondence_distance`.
111///    c. Build the weighted 6×6 Gauss-Newton system using combined covariances
112///       `M_i = C_i^T + R C_i^S R^T`.
113///    d. Solve for the incremental 6-DOF update and compose with the current transform.
114///    e. Declare convergence when |ΔMSE| < `convergence_threshold`.
115///
116/// # References
117/// Segal, A., Haehnel, D., & Thrun, S. (2009). *Generalized-ICP.*
118/// Robotics: Science and Systems.
119pub fn gicp(
120    source: &PointCloud<Point3f>,
121    target: &PointCloud<Point3f>,
122    init: Isometry3<f32>,
123    config: GicpConfig,
124) -> Result<ICPResult> {
125    if source.is_empty() || target.is_empty() {
126        return Err(Error::InvalidData(
127            "GICP: source or target point cloud is empty".into(),
128        ));
129    }
130    if config.max_iterations == 0 {
131        return Err(Error::InvalidData(
132            "GICP: max_iterations must be > 0".into(),
133        ));
134    }
135    // Covariance estimation needs at least k neighbours per point; fewer points
136    // than k_correspondences produces a degenerate (global) covariance.
137    let min_k = config.k_correspondences.max(4);
138    if source.len() < min_k || target.len() < min_k {
139        return Err(Error::InvalidData(format!(
140            "GICP: clouds must have at least {} points for reliable covariance estimation \
141             (k_correspondences={}); got source={}, target={}",
142            min_k, config.k_correspondences, source.len(), target.len()
143        )));
144    }
145    // Reject degenerate clouds where all points lie on a plane or line — GICP
146    // requires 3-D structure to build meaningful per-point covariance matrices.
147    for (label, cloud) in [("source", source), ("target", target)] {
148        let mut min = [f32::INFINITY; 3];
149        let mut max = [f32::NEG_INFINITY; 3];
150        for p in &cloud.points {
151            for ax in 0..3 {
152                min[ax] = min[ax].min(p.coords[ax]);
153                max[ax] = max[ax].max(p.coords[ax]);
154            }
155        }
156        let min_extent = (0..3).map(|ax| max[ax] - min[ax]).fold(f32::INFINITY, f32::min);
157        if min_extent < 1e-4 {
158            return Err(Error::InvalidData(format!(
159                "GICP: {label} point cloud appears to be coplanar or collinear \
160                 (smallest bounding-box dimension = {min_extent:.2e}); \
161                 GICP requires 3-D structure"
162            )));
163        }
164    }
165
166    // Step 1: per-point covariances.
167    let source_covs = compute_covariances(&source.points, config.k_correspondences)?;
168    let target_covs = compute_covariances(&target.points, config.k_correspondences)?;
169
170    // Step 2: KD-tree over target.
171    let target_tree = KdTree::new(&target.points)?;
172
173    let mut current_transform = init;
174    let mut prev_mse = f32::INFINITY;
175    let mut final_corr: Vec<(usize, usize)> = Vec::new();
176
177    for iteration in 0..config.max_iterations {
178        // Transform source points with the current estimate.
179        let transformed_source: Vec<Point3f> =
180            source.points.iter().map(|p| current_transform * p).collect();
181
182        // Rotation matrix for covariance transformation.
183        let r_mat = *current_transform.rotation.to_rotation_matrix().matrix();
184
185        // Accumulate the 6×6 Gauss-Newton Hessian H and gradient g.
186        // Pose parameterisation: x = [δω_x, δω_y, δω_z, δt_x, δt_y, δt_z]
187        // (rotation first, matching the existing point-to-plane convention).
188        let mut h = Matrix6::<f32>::zeros();
189        let mut gvec = Vector6::<f32>::zeros();
190        let mut n_corr = 0usize;
191        let mut mse_sum = 0.0f32;
192        let mut corr_pairs: Vec<(usize, usize)> = Vec::new();
193
194        for (src_idx, (ts, c_s)) in
195            transformed_source.iter().zip(source_covs.iter()).enumerate()
196        {
197            let neighbours = target_tree.find_k_nearest(ts, 1);
198            if neighbours.is_empty() {
199                continue;
200            }
201            let (tgt_idx, dist) = neighbours[0];
202            if dist > config.max_correspondence_distance {
203                continue;
204            }
205
206            let c_t = &target_covs[tgt_idx];
207
208            // Combined covariance M = C_t + R C_s R^T.
209            let m = c_t + r_mat * c_s * r_mat.transpose();
210            let m_inv = match m.try_inverse() {
211                Some(inv) => inv,
212                None => continue, // numerically degenerate — skip
213            };
214
215            // Residual r = t_i – T s_i.
216            let residual = target.points[tgt_idx].coords - ts.coords;
217
218            // 3×6 Jacobian J = [A | I] where A = −skew(T s_i).
219            // Contribution to H and g:
220            //   H_rr (0..3, 0..3) = A^T M⁻¹ A
221            //   H_rt (0..3, 3..6) = A^T M⁻¹
222            //   H_tr (3..6, 0..3) = M⁻¹ A = H_rt^T
223            //   H_tt (3..6, 3..6) = M⁻¹
224            //   g_r  (0..3)       = A^T M⁻¹ r
225            //   g_t  (3..6)       = M⁻¹ r
226
227            let a = -skew_sym(&ts.coords); // rotation block
228            let h_rr = a.transpose() * (m_inv * a);
229            let h_rt = a.transpose() * m_inv;
230            let wr = m_inv * residual;
231            let g_r = a.transpose() * wr;
232
233            for i in 0..3 {
234                for j in 0..3 {
235                    h[(i, j)] += h_rr[(i, j)];
236                    h[(i, j + 3)] += h_rt[(i, j)];
237                    h[(i + 3, j)] += h_rt[(j, i)]; // H_tr = H_rt^T
238                    h[(i + 3, j + 3)] += m_inv[(i, j)];
239                }
240                gvec[i] += g_r[i];
241                gvec[i + 3] += wr[i];
242            }
243
244            n_corr += 1;
245            mse_sum += dist * dist;
246            corr_pairs.push((src_idx, tgt_idx));
247        }
248
249        if n_corr < 6 {
250            return Err(Error::Algorithm(
251                "GICP: insufficient correspondences (need ≥ 6)".into(),
252            ));
253        }
254
255        let mse = mse_sum / n_corr as f32;
256
257        // Solve H δ = g (Cholesky preferred; LU fallback).
258        let delta = if let Some(chol) = h.cholesky() {
259            chol.solve(&gvec)
260        } else {
261            h.lu().solve(&gvec).ok_or_else(|| {
262                Error::Algorithm("GICP: Gauss-Newton system is ill-conditioned".into())
263            })?
264        };
265
266        // Compose incremental transform δT from x = [δω_x, δω_y, δω_z, δt_x, δt_y, δt_z].
267        let delta_rot = UnitQuaternion::from_axis_angle(&Vector3f::z_axis(), delta[2])
268            * UnitQuaternion::from_axis_angle(&Vector3f::y_axis(), delta[1])
269            * UnitQuaternion::from_axis_angle(&Vector3f::x_axis(), delta[0]);
270        let delta_iso = Isometry3::from_parts(
271            Translation3::new(delta[3], delta[4], delta[5]),
272            delta_rot,
273        );
274        current_transform = delta_iso * current_transform;
275
276        // Convergence check.
277        if (prev_mse - mse).abs() < config.convergence_threshold {
278            return Ok(ICPResult {
279                transformation: current_transform,
280                mse,
281                iterations: iteration + 1,
282                converged: true,
283                correspondences: corr_pairs,
284            });
285        }
286
287        prev_mse = mse;
288        final_corr = corr_pairs;
289    }
290
291    Ok(ICPResult {
292        transformation: current_transform,
293        mse: prev_mse,
294        iterations: config.max_iterations,
295        converged: false,
296        correspondences: final_corr,
297    })
298}
299
300// ---------------------------------------------------------------------------
301// Tests
302// ---------------------------------------------------------------------------
303
304#[cfg(test)]
305mod tests {
306    use super::*;
307    use approx::assert_relative_eq;
308    #[allow(unused_imports)]
309    use nalgebra::Translation3;
310
311    /// Fibonacci-sphere point cloud with outward unit normals (not used by GICP
312    /// directly, but a good well-conditioned test surface).
313    fn make_sphere(n: usize, radius: f32) -> PointCloud<Point3f> {
314        let mut cloud = PointCloud::new();
315        let golden = std::f32::consts::PI * (3.0 - 5.0_f32.sqrt());
316        for i in 0..n {
317            let y = 1.0 - (i as f32 / (n as f32 - 1.0).max(1.0)) * 2.0;
318            let r = (1.0 - y * y).max(0.0_f32).sqrt();
319            let theta = golden * i as f32;
320            cloud.push(Point3f::new(theta.cos() * r * radius, y * radius, theta.sin() * r * radius));
321        }
322        cloud
323    }
324
325    #[test]
326    fn gicp_identity_converges() {
327        let cloud = make_sphere(100, 3.0);
328        let config = GicpConfig { max_iterations: 30, ..Default::default() };
329        let result = gicp(&cloud, &cloud, Isometry3::identity(), config).unwrap();
330        assert!(result.converged, "should converge for identical clouds");
331        assert!(result.mse < 1e-4, "mse={}", result.mse);
332    }
333
334    #[test]
335    fn gicp_recovers_small_translation() {
336        let source = make_sphere(150, 3.0);
337        let shift = Vector3f::new(0.1, 0.0, 0.0);
338        let target = PointCloud::from_points(source.points.iter().map(|p| p + shift).collect());
339
340        let config = GicpConfig {
341            max_iterations: 60,
342            max_correspondence_distance: 2.0,
343            ..Default::default()
344        };
345        let result = gicp(&source, &target, Isometry3::identity(), config).unwrap();
346
347        let t_err = (result.transformation.translation.vector - shift).magnitude();
348        assert!(t_err < 0.05, "translation error={}", t_err);
349        assert!(result.mse < 0.1, "mse={}", result.mse);
350    }
351
352    #[test]
353    fn gicp_empty_source_errors() {
354        let empty: PointCloud<Point3f> = PointCloud::new();
355        let cloud = make_sphere(30, 1.0);
356        assert!(gicp(&empty, &cloud, Isometry3::identity(), GicpConfig::default()).is_err());
357    }
358
359    #[test]
360    fn gicp_zero_iterations_errors() {
361        let cloud = make_sphere(30, 1.0);
362        let config = GicpConfig { max_iterations: 0, ..Default::default() };
363        assert!(gicp(&cloud, &cloud, Isometry3::identity(), config).is_err());
364    }
365
366    #[test]
367    fn gicp_result_fields_populated() {
368        let cloud = make_sphere(60, 2.0);
369        let config = GicpConfig { max_iterations: 10, ..Default::default() };
370        let result = gicp(&cloud, &cloud, Isometry3::identity(), config).unwrap();
371        assert!(result.iterations > 0);
372        assert!(!result.correspondences.is_empty());
373    }
374
375    // -------------------------------------------------------------------------
376    // Rotation recovery
377    //
378    // ICP is a LOCAL optimizer: it refines an initial guess rather than solving
379    // the global alignment problem.  Tests use either:
380    //   (a) identity init + tiny rotation (< half the average point spacing), or
381    //   (b) a near-correct initial guess for larger rotations, matching the real
382    //       SLAM use-case where IMU/odometry supplies a rough estimate.
383    // -------------------------------------------------------------------------
384
385    #[test]
386    fn gicp_recovers_tiny_rotation_from_identity() {
387        // 300-point sphere → average angular spacing ≈ 11°.
388        // 2° rotation chord ≈ 0.21 m — well within the half-spacing convergence
389        // basin so every nearest-neighbour correspondence is the correct one.
390        let source = make_sphere(300, 3.0);
391        let angle = 2_f32.to_radians();
392        let rot = UnitQuaternion::from_axis_angle(&Vector3f::z_axis(), angle);
393        let target = PointCloud::from_points(
394            source.points.iter().map(|p| rot * p).collect(),
395        );
396
397        let config = GicpConfig {
398            max_iterations: 60,
399            max_correspondence_distance: 0.8,
400            ..Default::default()
401        };
402        let result = gicp(&source, &target, Isometry3::identity(), config).unwrap();
403
404        let rot_err = result.transformation.rotation.angle_to(&rot);
405        assert!(
406            rot_err < 0.5_f32.to_radians(),
407            "rotation error = {:.2}°",
408            rot_err.to_degrees()
409        );
410        assert!(result.mse < 0.01, "mse={}", result.mse);
411    }
412
413    #[test]
414    fn gicp_refines_rotation_from_near_correct_init() {
415        // Simulate the SLAM use-case: odometry gives an 80%-correct initial guess
416        // (6° out of 8°).  GICP refines the 2° residual to < 0.5°.
417        let source = make_sphere(200, 3.0);
418        let true_angle = 8_f32.to_radians();
419        let init_angle = 6_f32.to_radians(); // ← provided by odometry / IMU
420        let true_rot = UnitQuaternion::from_axis_angle(&Vector3f::z_axis(), true_angle);
421        let target = PointCloud::from_points(
422            source.points.iter().map(|p| true_rot * p).collect(),
423        );
424        let init = Isometry3::from_parts(
425            Translation3::identity(),
426            UnitQuaternion::from_axis_angle(&Vector3f::z_axis(), init_angle),
427        );
428
429        let config = GicpConfig {
430            max_iterations: 60,
431            max_correspondence_distance: 0.8,
432            ..Default::default()
433        };
434        let result = gicp(&source, &target, init, config).unwrap();
435
436        let rot_err = result.transformation.rotation.angle_to(&true_rot);
437        assert!(
438            rot_err < 0.5_f32.to_radians(),
439            "rotation error = {:.2}°",
440            rot_err.to_degrees()
441        );
442    }
443
444    #[test]
445    fn gicp_refines_combined_rotation_and_translation() {
446        // Near-correct init for combined 6° rotation + 0.3 m translation.
447        let source = make_sphere(200, 3.0);
448        let true_angle = 6_f32.to_radians();
449        let true_shift = Vector3f::new(0.3, 0.0, 0.0);
450        let true_rot = UnitQuaternion::from_axis_angle(&Vector3f::y_axis(), true_angle);
451        let true_iso = Isometry3::from_parts(
452            Translation3::new(true_shift.x, true_shift.y, true_shift.z),
453            true_rot,
454        );
455        let target = PointCloud::from_points(
456            source.points.iter().map(|p| true_iso * p).collect(),
457        );
458
459        // Initial guess: 80% of the rotation, 80% of the translation.
460        let init = Isometry3::from_parts(
461            Translation3::new(0.24, 0.0, 0.0),
462            UnitQuaternion::from_axis_angle(&Vector3f::y_axis(), 4.8_f32.to_radians()),
463        );
464
465        let config = GicpConfig {
466            max_iterations: 80,
467            max_correspondence_distance: 0.8,
468            ..Default::default()
469        };
470        let result = gicp(&source, &target, init, config).unwrap();
471
472        let t_err = (result.transformation.translation.vector - true_shift).magnitude();
473        let rot_err = result.transformation.rotation.angle_to(&true_rot);
474        assert!(t_err < 0.05, "translation error={}", t_err);
475        assert!(
476            rot_err < 0.5_f32.to_radians(),
477            "rotation error = {:.2}°",
478            rot_err.to_degrees()
479        );
480    }
481
482    // -------------------------------------------------------------------------
483    // Noise robustness
484    // -------------------------------------------------------------------------
485
486    #[test]
487    fn gicp_robust_to_gaussian_noise() {
488        // Target = source + uniform ±0.05 m noise per axis.
489        let source = make_sphere(200, 3.0);
490        let noise_half = 0.05_f32;
491        // Deterministic noise via a simple LCG so tests are reproducible.
492        let target = PointCloud::from_points(
493            source.points.iter().enumerate().map(|(i, p)| {
494                let t = (i as f32 * 1.6180339887_f32).sin() * noise_half;
495                let u = (i as f32 * 2.7182818284_f32).cos() * noise_half;
496                let v = (i as f32 * 3.1415926535_f32).sin() * noise_half;
497                Point3f::new(p.x + t, p.y + u, p.z + v)
498            }).collect(),
499        );
500
501        let config = GicpConfig {
502            max_iterations: 50,
503            max_correspondence_distance: 1.0,
504            ..Default::default()
505        };
506        let result = gicp(&source, &target, Isometry3::identity(), config).unwrap();
507
508        // MSE should be in the noise-variance ballpark (3 × 0.05² ≈ 0.0075),
509        // well below the correspondence-distance threshold.
510        assert!(result.mse < 0.05, "mse={}", result.mse);
511        let t_err = result.transformation.translation.vector.magnitude();
512        assert!(t_err < 0.1, "spurious translation drift={}", t_err);
513    }
514
515    #[test]
516    fn gicp_robust_to_outlier_points() {
517        // Target = source union 20 far outliers. Tight correspondence distance
518        // should exclude the outliers so the result stays close to identity.
519        let source = make_sphere(200, 3.0);
520        let mut target_pts = source.points.clone();
521        for i in 0..20 {
522            let t = i as f32;
523            target_pts.push(Point3f::new(t * 7.3 - 50.0, t * 3.1 - 30.0, t * 5.7 - 40.0));
524        }
525        let target = PointCloud::from_points(target_pts);
526
527        let config = GicpConfig {
528            max_iterations: 40,
529            max_correspondence_distance: 0.5,
530            ..Default::default()
531        };
532        let result = gicp(&source, &target, Isometry3::identity(), config).unwrap();
533
534        let t_err = result.transformation.translation.vector.magnitude();
535        assert!(t_err < 0.05, "spurious drift from outliers={}", t_err);
536        assert!(result.mse < 0.01);
537    }
538
539    // -------------------------------------------------------------------------
540    // Degenerate inputs
541    // -------------------------------------------------------------------------
542
543    #[test]
544    fn gicp_too_few_points_errors() {
545        // Fewer points than k_correspondences (default 20).
546        let cloud = make_sphere(10, 1.0);
547        assert!(
548            gicp(&cloud, &cloud, Isometry3::identity(), GicpConfig::default()).is_err(),
549            "expected error for cloud with too few points"
550        );
551    }
552
553    #[test]
554    fn gicp_coplanar_points_errors() {
555        // All points on Z = 0 plane — bounding box extent in Z is 0.
556        let pts: Vec<Point3f> = (0..50)
557            .flat_map(|i| (0..50).map(move |j| Point3f::new(i as f32 * 0.1, j as f32 * 0.1, 0.0)))
558            .collect();
559        let cloud = PointCloud::from_points(pts);
560        assert!(
561            gicp(&cloud, &cloud, Isometry3::identity(), GicpConfig::default()).is_err(),
562            "expected error for coplanar point cloud"
563        );
564    }
565}