mapping_algorithms/point_clouds/icp/
mod.rs

1// SPDX-License-Identifier: MIT
2/*
3 * Copyright (c) [2023 - Present] Emily Matheys <emilymatt96@gmail.com>
4 *
5 * Permission is hereby granted, free of charge, to any person obtaining a copy
6 * of this software and associated documentation files (the "Software"), to deal
7 * in the Software without restriction, including without limitation the rights
8 * to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
9 * copies of the Software, and to permit persons to whom the Software is
10 * furnished to do so, subject to the following conditions:
11 *
12 * The above copyright notice and this permission notice shall be included in all
13 * copies or substantial portions of the Software.
14 *
15 * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
16 * IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
17 * FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
18 * AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
19 * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
20 * OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
21 * SOFTWARE.
22 */
23
24pub use types::{ICPConfiguration, ICPConfigurationBuilder, ICPError, ICPResult, ICPSuccess};
25
26use nalgebra::{ComplexField, Isometry, Point, RealField, SimdRealField};
27use num_traits::{AsPrimitive, Bounded};
28
29use crate::{
30    kd_tree::KDTree,
31    point_clouds::find_nearest_neighbour_naive,
32    types::{AbstractIsometry, IsNan, IsometryAbstractor},
33    Sum, Vec,
34};
35
36use helpers::{calculate_mse, get_rotation_matrix_and_centroids};
37
38mod helpers;
39mod types;
40
41/// A single iteration of the ICP function, allowing for any input and output, usually used for debugging or visualization
42///
43/// # Arguments
44/// * `points_a`: A slice of [`Point`], representing the source point cloud.
45/// * `transformed_points`: A mutable slice of [`Point`], representing the transformed source point cloud, this will be transformed further by the function.
46/// * `points_b`: A slice of [`Point`], representing the target point cloud.
47/// * `target_points_tree`: An [`Option<KDTree<T, N>>`], this is usually created by the ICP function if `config.use_kd` is `true`
48/// * `current_transform`: A mutable reference to the [`Isometry`] used to transform the source points, this will gradually change with each iteration.
49/// * `current_mse`: A mutable reference of a `T`, this will be updated by the function to the latest MSE, which is then used by the ICP function to determine an exit strategy.
50/// * `config`: a reference to an [`ICPConfiguration`], specifying the behaviour of the algorithm.
51///
52/// # Generics
53/// * `T`: Either [`prim@f32`] or [`prim@f64`].
54/// * `R`: Either a [`UnitComplex`](nalgebra::UnitComplex) or a [`UnitQuaternion`](nalgebra::UnitQuaternion) of `T`, depednding on `N`.
55/// * `N`: a usize, either `2` or `3`.
56///
57/// # Returns
58/// An [`ICPSuccess`] struct with an [`Isometry`] transform with a `T` precision, or an error message explaining what went wrong.
59///
60/// [^convergence_note]: This does not guarantee that the transformation is correct, only that no further benefit can be gained by running another iteration.
61#[cfg_attr(
62    feature = "tracing",
63    tracing::instrument("ICP Algorithm Iteration", skip_all, level = "info")
64)]
65pub fn icp_iteration<T, const N: usize>(
66    points_a: &[Point<T, N>],
67    transformed_points: &mut [Point<T, N>],
68    points_b: &[Point<T, N>],
69    target_points_tree: Option<&KDTree<T, N>>,
70    current_transform: &mut Isometry<
71        T,
72        <IsometryAbstractor<T, N> as AbstractIsometry<T, N>>::RotType,
73        N,
74    >,
75    current_mse: &mut T,
76    config: &ICPConfiguration<T>,
77) -> Result<T, ICPError<T, N>>
78where
79    T: Bounded + Copy + Default + RealField + Sum + SimdRealField,
80    usize: AsPrimitive<T>,
81    IsometryAbstractor<T, N>: AbstractIsometry<T, N>,
82{
83    let closest_points = transformed_points.iter().try_fold(
84        Vec::with_capacity(transformed_points.len()),
85        |mut accumulator, transformed_point_a| {
86            accumulator.push(
87                target_points_tree
88                    .and_then(|kd_tree| kd_tree.nearest(transformed_point_a))
89                    .or_else(|| find_nearest_neighbour_naive(transformed_point_a, points_b))
90                    .ok_or(ICPError::NoNearestNeighbour)?,
91            );
92
93            Ok(accumulator)
94        },
95    )?;
96
97    let (rot_mat, mean_a, mean_b) =
98        get_rotation_matrix_and_centroids(transformed_points, &closest_points);
99
100    *current_transform =
101        IsometryAbstractor::<T, N>::update_transform(current_transform, mean_a, mean_b, &rot_mat);
102
103    for (idx, point_a) in points_a.iter().enumerate() {
104        transformed_points[idx] = current_transform.transform_point(point_a);
105    }
106    let new_mse = calculate_mse(transformed_points, closest_points.as_slice());
107    log::trace!("New MSE: {new_mse}");
108
109    // If the MSE difference is lower than the threshold, then this is as good as it gets
110    if config
111        .mse_absolute_threshold
112        .map(|thres| new_mse < thres)
113        .unwrap_or_default()
114        || <T as ComplexField>::abs(*current_mse - new_mse) < config.mse_interval_threshold
115    {
116        return Ok(new_mse);
117    }
118
119    *current_mse = new_mse;
120    Err(ICPError::IterationDidNotConverge((mean_a, mean_b)))
121}
122
123/// A free-form version of the ICP function, allowing for any input and output, under the constraints of the function
124///
125/// # Arguments
126/// * `points_a`: A slice of [`Point`], representing the source point cloud.
127/// * `points_b`: A slice of [`Point`], representing the target point cloud.
128/// * `config`: a reference to an [`ICPConfiguration`], specifying the behaviour of the algorithm.
129///
130/// # Generics
131/// * `T`: Either [`prim@f32`] or [`prim@f64`]
132/// * `R`: An ['SVDIsometry'] of `T` and `N`
133/// * `N`: a usize, either `2` or `3`
134///
135/// # Returns
136/// An [`ICPSuccess`] struct with an [`Isometry`] transform with a `T` precision, or an error message explaining what went wrong.
137///
138/// [^convergence_note]: This does not guarantee that the transformation is correct, only that no further benefit can be gained by running another iteration.
139#[cfg_attr(
140    feature = "tracing",
141    tracing::instrument("Full ICP Algorithm", skip_all, level = "info")
142)]
143pub fn icp<T, const N: usize>(
144    points_a: &[Point<T, N>],
145    points_b: &[Point<T, N>],
146    config: ICPConfiguration<T>,
147) -> ICPResult<T, <IsometryAbstractor<T, N> as AbstractIsometry<T, N>>::RotType, N>
148where
149    T: Bounded + Copy + Default + IsNan + RealField + Sum,
150    usize: AsPrimitive<T>,
151    IsometryAbstractor<T, N>: AbstractIsometry<T, N>,
152{
153    if points_a.is_empty() {
154        return Err(ICPError::SourcePointCloudEmpty);
155    }
156
157    if points_b.is_empty() {
158        return Err(ICPError::TargetPointCloudEmpty);
159    }
160
161    if config.max_iterations == 0 {
162        return Err(ICPError::IterationNumIsZero);
163    }
164
165    if config.mse_interval_threshold <= T::default_epsilon() {
166        return Err(ICPError::MSEIntervalThreshold);
167    }
168
169    if config
170        .mse_absolute_threshold
171        .map(|thres| thres.is_nan() || thres <= T::default_epsilon())
172        .unwrap_or_default()
173    {
174        return Err(ICPError::MSEAbsoluteThreshold);
175    }
176
177    let mut points_to_transform = points_a.to_vec();
178    let target_points_tree = config.use_kd_tree.then_some(KDTree::from(points_b));
179    let mut current_transform = Isometry::identity();
180    let mut current_mse = <T as Bounded>::max_value();
181
182    for iteration_num in 0..config.max_iterations {
183        log::trace!(
184            "Running iteration number {iteration_num}/{}",
185            config.max_iterations
186        );
187        if let Ok(mse) = icp_iteration::<T, N>(
188            points_a,
189            &mut points_to_transform,
190            points_b,
191            target_points_tree.as_ref(),
192            &mut current_transform,
193            &mut current_mse,
194            &config,
195        ) {
196            log::trace!("Converged after {iteration_num} iterations with an MSE of {mse}");
197            return Ok(ICPSuccess {
198                transform: current_transform,
199                mse,
200                iteration_num,
201            });
202        }
203    }
204
205    Err(ICPError::AlrogithmDidNotConverge)
206}
207
208#[cfg(feature = "pregenerated")]
209macro_rules! impl_icp_algorithm {
210    ($precision:expr, $doc:tt, $nd:expr, $rot_type:expr) => {
211        ::paste::paste! {
212            #[doc = "A premade variant of the ICP algorithm function, in " $nd "D space and " $doc "-precision floats."]
213            pub fn [<icp_$nd d>](points_a: &[Point<$precision, $nd>],
214                points_b: &[Point<$precision, $nd>],
215                config: ICPConfiguration<$precision>) -> ICPResult<$precision, $rot_type<$precision>, $nd> {
216                    super::icp(points_a, points_b, config)
217            }
218        }
219    };
220
221    ($precision:expr, doc $doc:tt) => {
222        ::paste::paste! {
223            pub(super) mod [<$doc _precision>] {
224                use nalgebra::{Point, UnitComplex, UnitQuaternion};
225                use super::{ICPConfiguration, ICPResult};
226
227                impl_icp_algorithm!($precision, $doc, 2, UnitComplex);
228                impl_icp_algorithm!($precision, $doc, 3, UnitQuaternion);
229            }
230        }
231    }
232}
233
234#[cfg(feature = "pregenerated")]
235impl_icp_algorithm!(f32, doc single);
236#[cfg(feature = "pregenerated")]
237impl_icp_algorithm!(f64, doc double);
238
239#[cfg(test)]
240mod tests {
241    use nalgebra::{Isometry2, Isometry3, UnitComplex, Vector2, Vector3};
242
243    use crate::{
244        array,
245        point_clouds::{generate_point_cloud, transform_point_cloud},
246    };
247
248    use super::*;
249
250    #[test]
251    fn test_icp_errors() {
252        let points = generate_point_cloud(10, array::from_fn(|_| -15.0..=15.0));
253        let config_builder = ICPConfiguration::builder();
254
255        let mut res: ICPResult<f32, UnitComplex<f32>, 2> =
256            icp(&[], points.as_slice(), config_builder.build());
257        assert_eq!(res.unwrap_err(), ICPError::SourcePointCloudEmpty);
258
259        res = icp(points.as_slice(), &[], config_builder.build());
260        assert_eq!(res.unwrap_err(), ICPError::TargetPointCloudEmpty);
261
262        res = icp(
263            points.as_slice(),
264            points.as_slice(),
265            config_builder.with_max_iterations(0).build(),
266        );
267        assert_eq!(res.unwrap_err(), ICPError::IterationNumIsZero);
268
269        res = icp(
270            points.as_slice(),
271            points.as_slice(),
272            config_builder.with_mse_interval_threshold(0.0).build(),
273        );
274        assert_eq!(res.unwrap_err(), ICPError::MSEIntervalThreshold);
275
276        res = icp(
277            points.as_slice(),
278            points.as_slice(),
279            config_builder
280                .with_absolute_mse_threshold(Some(0.0))
281                .build(),
282        );
283        assert_eq!(res.unwrap_err(), ICPError::MSEAbsoluteThreshold);
284    }
285
286    #[test]
287    fn test_no_convegence() {
288        let points = generate_point_cloud(1000, array::from_fn(|_| -15.0..=15.0));
289        let translation = Vector2::new(-12.5, 7.3);
290        let isom = Isometry2::new(translation, 90.0f32.to_radians());
291        let points_transformed = transform_point_cloud(&points, isom);
292
293        let res = icp(
294            points.as_slice(),
295            points_transformed.as_slice(),
296            ICPConfiguration::builder()
297                .with_max_iterations(1) // No chance something like this could converge, and definitely not in 1 iteration
298                .with_mse_interval_threshold(0.001)
299                .build(),
300        );
301        assert!(res.is_err());
302        assert_eq!(res.unwrap_err(), ICPError::AlrogithmDidNotConverge);
303    }
304
305    #[test]
306    // This is for code coverage purposes, ensure that absolute MSE is used instead of interval
307    fn test_icp_absolute_threshold() {
308        let points = generate_point_cloud(100, array::from_fn(|_| -15.0..=15.0));
309        let translation = Vector2::new(-0.8, 1.3);
310        let isom = Isometry2::new(translation, 0.1);
311        let points_transformed = transform_point_cloud(&points, isom);
312
313        let res = icp(
314            points.as_slice(),
315            points_transformed.as_slice(),
316            ICPConfiguration::builder()
317                .with_max_iterations(10)
318                .with_absolute_mse_threshold(Some(0.1))
319                .with_mse_interval_threshold(0.001)
320                .build(),
321        );
322        assert!(res.is_ok());
323        assert!(res.unwrap().mse < 0.1);
324    }
325
326    #[test]
327    fn test_icp_2d() {
328        let points = generate_point_cloud(100, array::from_fn(|_| -15.0..=15.0));
329        let translation = Vector2::new(-0.8, 1.3);
330        let isom = Isometry2::new(translation, 0.1);
331        let points_transformed = transform_point_cloud(&points, isom);
332
333        let res = icp(
334            points.as_slice(),
335            points_transformed.as_slice(),
336            ICPConfiguration::builder()
337                .with_max_iterations(10)
338                .with_mse_interval_threshold(0.01)
339                .build(),
340        );
341        assert!(res.is_ok());
342        assert!(res.unwrap().mse < 0.01);
343    }
344
345    #[test]
346    fn test_icp_2d_with_kd() {
347        let points = generate_point_cloud(100, array::from_fn(|_| -15.0..=15.0));
348        let isom = Isometry2::new(Vector2::new(-0.8, 1.3), 0.1);
349        let points_transformed = transform_point_cloud(&points, isom);
350
351        let res = icp(
352            points.as_slice(),
353            points_transformed.as_slice(),
354            ICPConfiguration::builder()
355                .with_kd_tree(true)
356                .with_max_iterations(50)
357                .with_mse_interval_threshold(0.01)
358                .build(),
359        );
360        assert!(res.is_ok());
361        assert!(res.unwrap().mse < 0.01);
362    }
363
364    #[test]
365    fn test_icp_3d() {
366        let points = generate_point_cloud(500, array::from_fn(|_| -15.0..=15.0));
367        let translation = Vector3::new(-0.8, 1.3, 0.2);
368        let rotation = Vector3::new(0.1, 0.2, -0.21);
369        let isom = Isometry3::new(translation, rotation);
370        let points_transformed = transform_point_cloud(&points, isom);
371
372        let res = icp(
373            points.as_slice(),
374            points_transformed.as_slice(),
375            ICPConfiguration::builder()
376                .with_max_iterations(50)
377                .with_mse_interval_threshold(0.01)
378                .build(),
379        );
380        assert!(res.is_ok());
381        assert!(res.unwrap().mse < 0.05);
382    }
383
384    #[test]
385    fn test_icp_3d_with_kd() {
386        let points = generate_point_cloud(500, array::from_fn(|_| -15.0..=15.0));
387        let translation = Vector3::new(-0.8, 1.3, 0.2);
388        let rotation = Vector3::new(0.1, 0.2, -0.21);
389        let isom = Isometry3::new(translation, rotation);
390        let points_transformed = transform_point_cloud(&points, isom);
391
392        let res = icp(
393            points.as_slice(),
394            points_transformed.as_slice(),
395            ICPConfiguration::builder()
396                .with_kd_tree(true)
397                .with_max_iterations(50)
398                .with_mse_interval_threshold(0.01)
399                .build(),
400        );
401        assert!(res.is_ok());
402        assert!(res.unwrap().mse < 0.05);
403    }
404}