Skip to main content

hoomd_geometry/shape/
convex_polytope.rs

1// Copyright (c) 2024-2026 The Regents of the University of Michigan.
2// Part of hoomd-rs, released under the BSD 3-Clause License.
3
4//! N-Dimensional generalization of a convex polyhedron.
5
6use serde::{Deserialize, Serialize};
7
8use crate::{BoundingSphereRadius, Error, SupportMapping};
9use hoomd_utility::valid::PositiveReal;
10use hoomd_vector::{Cartesian, InnerProduct};
11
12/// A faceted solid defined by the convex hull of its vertices.
13///
14/// # Examples
15///
16/// Construction and basic methods:
17/// ```
18/// use approxim::assert_relative_eq;
19/// use hoomd_geometry::{BoundingSphereRadius, shape::ConvexPolyhedron};
20///
21/// # fn main() -> Result<(), hoomd_geometry::Error> {
22/// let tetrahedron = ConvexPolyhedron::with_vertices(vec![
23///     [1.0, 1.0, 1.0].into(),
24///     [1.0, -1.0, -1.0].into(),
25///     [-1.0, 1.0, -1.0].into(),
26///     [-1.0, -1.0, 1.0].into(),
27/// ])?;
28///
29/// let bounding_radius = tetrahedron.bounding_sphere_radius();
30///
31/// assert_relative_eq!(bounding_radius.get(), 3.0_f64.sqrt());
32/// # Ok(())
33/// # }
34/// ```
35///
36/// Intersection tests:
37/// ```
38/// use hoomd_geometry::{Convex, IntersectsAt, shape::ConvexPolygon};
39/// use hoomd_vector::{Angle, Cartesian};
40/// use std::f64::consts::PI;
41///
42/// # fn main() -> Result<(), hoomd_geometry::Error> {
43/// let rectangle = ConvexPolygon::with_vertices([
44///     [-2.0, -1.0].into(),
45///     [2.0, -1.0].into(),
46///     [2.0, 1.0].into(),
47///     [-2.0, 1.0].into(),
48/// ])?;
49/// let rectangle = Convex(rectangle);
50///
51/// assert!(!rectangle.intersects_at(
52///     &rectangle,
53///     &[0.0, 2.1].into(),
54///     &Angle::default()
55/// ));
56/// assert!(rectangle.intersects_at(
57///     &rectangle,
58///     &[0.0, 2.1].into(),
59///     &Angle::from(PI / 2.0)
60/// ));
61/// # Ok(())
62/// # }
63/// ```
64#[derive(Clone, Debug, PartialEq, Serialize, Deserialize)]
65pub struct ConvexPolytope<const N: usize> {
66    /// The vertices of the shape.
67    vertices: Vec<Cartesian<N>>,
68    /// The radius of a bounding sphere of the geometry.
69    pub(crate) bounding_radius: PositiveReal,
70}
71
72/// A faceted convex body in two dimensions.
73///
74/// ```rust
75/// use hoomd_geometry::shape::ConvexPolygon;
76///
77/// # fn main() -> Result<(), hoomd_geometry::Error> {
78/// let hexagon = ConvexPolygon::regular(6);
79/// let square = ConvexPolygon::with_vertices([
80///     [-1.0, -1.0].into(),
81///     [1.0, -1.0].into(),
82///     [1.0, 1.0].into(),
83///     [-1.0, 1.0].into(),
84/// ])?;
85/// # Ok(())
86/// # }
87/// ```
88pub type ConvexPolygon = ConvexPolytope<2>;
89
90/// A faceted convex body in three dimensions.
91///
92/// # Example
93///
94/// ```
95/// use hoomd_geometry::shape::{ConvexPolyhedron, Simplex3};
96/// # fn main() -> Result<(), hoomd_geometry::Error> {
97/// let poly = ConvexPolyhedron::with_vertices(vec![
98///     [1.0, 1.0, 1.0].into(),
99///     [1.0, -1.0, -1.0].into(),
100///     [-1.0, 1.0, -1.0].into(),
101///     [-1.0, -1.0, 1.0].into(),
102/// ])?;
103///
104/// assert_eq!(poly.vertices(), Simplex3::default().vertices());
105/// # Ok(())
106/// # }
107/// ```
108pub type ConvexPolyhedron = ConvexPolytope<3>;
109
110impl ConvexPolytope<2> {
111    /// Create a regular *n*-gon with *n* vertices and circumradius 0.5.
112    ///
113    /// # Example
114    /// ```
115    /// use hoomd_geometry::shape::ConvexPolytope;
116    ///
117    /// let equilateral_triangle = ConvexPolytope::regular(3);
118    /// ```
119    #[inline]
120    #[must_use]
121    #[expect(
122        clippy::missing_panics_doc,
123        reason = "panic will never occur on a hard-coded constant"
124    )]
125    pub fn regular(n: usize) -> ConvexPolytope<2> {
126        ConvexPolytope {
127            vertices: (0..n)
128                .map(|x| {
129                    let theta = 2.0 * std::f64::consts::PI * (x as f64) / (n as f64);
130                    Cartesian::from([0.5 * f64::cos(theta), 0.5 * f64::sin(theta)])
131                })
132                .collect::<Vec<_>>(),
133            bounding_radius: 0.5
134                .try_into()
135                .expect("hard-coded constant should be positive"),
136        }
137    }
138}
139
140impl<const N: usize> ConvexPolytope<N> {
141    /// Create an `N`-polytope with the given vertices.
142    ///
143    /// # Example
144    /// ```
145    /// use hoomd_geometry::shape::ConvexPolytope;
146    ///
147    /// # fn main() -> Result<(), hoomd_geometry::Error> {
148    /// let equilateral_triangle = ConvexPolytope::with_vertices(vec![
149    ///     [1.0, 0.0].into(),
150    ///     [0.5, f64::sqrt(3.0) / 2.0].into(),
151    ///     [-0.5, f64::sqrt(3.0) / 2.0].into(),
152    /// ])?;
153    /// # Ok(())
154    /// # }
155    /// ```
156    /// # Errors
157    ///
158    /// `[Error::DegeneratePolytope]` when no vertices are provided.
159    #[inline]
160    pub fn with_vertices<I>(vertices: I) -> Result<ConvexPolytope<N>, Error>
161    where
162        I: IntoIterator<Item = Cartesian<N>>,
163    {
164        let vertices = vertices.into_iter().collect::<Vec<_>>();
165
166        if vertices.is_empty() {
167            return Err(Error::DegeneratePolytope);
168        }
169
170        Ok(ConvexPolytope {
171            bounding_radius: Self::bounding_radius(&vertices),
172            vertices,
173        })
174    }
175
176    /// The vertices of the shape.
177    #[inline]
178    #[must_use]
179    pub fn vertices(&self) -> &[Cartesian<N>] {
180        &self.vertices
181    }
182
183    /// Compute the bounding radius.
184    fn bounding_radius(vertices: &[Cartesian<N>]) -> PositiveReal {
185        vertices
186            .iter()
187            .map(Cartesian::norm_squared)
188            .fold(0.0, f64::max)
189            .sqrt()
190            .try_into()
191            .expect("convex polytope should have a positive bounding radius")
192    }
193}
194
195impl<const N: usize> SupportMapping<Cartesian<N>> for ConvexPolytope<N> {
196    #[inline]
197    fn support_mapping(&self, n: &Cartesian<N>) -> Cartesian<N> {
198        match N {
199            0 => Cartesian::<N>::default(),
200            1 => self.vertices[0],
201            _ => *self
202                .vertices
203                .iter()
204                .max_by(|a, b| {
205                    a.dot(n)
206                        .partial_cmp(&b.dot(n))
207                        .unwrap_or(std::cmp::Ordering::Equal)
208                })
209                .expect("the 0 match statement should handle empty vectors"),
210        }
211    }
212}
213
214impl<const N: usize> BoundingSphereRadius for ConvexPolytope<N> {
215    #[inline]
216    fn bounding_sphere_radius(&self) -> PositiveReal {
217        self.bounding_radius
218    }
219}
220
221#[cfg(test)]
222mod tests {
223    use super::*;
224    use crate::{Convex, IntersectsAt};
225    use hoomd_vector::{Angle, Cartesian, Rotate, Rotation, Versor};
226
227    use approxim::assert_relative_eq;
228    use rstest::*;
229    use std::f64::consts::{FRAC_1_SQRT_2, PI};
230
231    #[fixture]
232    fn simplex3() -> ConvexPolyhedron {
233        ConvexPolyhedron::with_vertices(vec![
234            [1.0, 1.0, 1.0].into(),
235            [1.0, -1.0, -1.0].into(),
236            [-1.0, 1.0, -1.0].into(),
237            [-1.0, -1.0, 1.0].into(),
238        ])
239        .unwrap()
240    }
241
242    #[fixture]
243    fn equilateral_triangle() -> ConvexPolytope<2> {
244        ConvexPolytope::with_vertices(vec![
245            [1.0, 0.0].into(),
246            [0.5, f64::sqrt(3.0) / 2.0].into(),
247            [-0.5, f64::sqrt(3.0) / 2.0].into(),
248        ])
249        .unwrap()
250    }
251
252    #[rstest]
253    fn test_bounding_radius_computed(
254        simplex3: ConvexPolyhedron,
255        equilateral_triangle: ConvexPolygon,
256    ) {
257        assert_eq!(simplex3.bounding_radius.get(), f64::sqrt(3.0));
258        assert_eq!(equilateral_triangle.bounding_radius.get(), f64::sqrt(1.0));
259    }
260
261    #[rstest]
262    fn test_bounding_radius_regular_polygons(#[values(1, 3, 8, 64)] n: usize) {
263        assert_eq!(ConvexPolygon::regular(n).bounding_radius.get(), 0.5);
264        assert_eq!(ConvexPolytope::regular(n).bounding_radius.get(), 0.5);
265    }
266
267    #[test]
268    fn degenerate_polytope() {
269        let result = ConvexPolytope::<3>::with_vertices([]);
270        assert_eq!(result, Err(Error::DegeneratePolytope));
271    }
272
273    #[test]
274    fn support_mapping_2d() {
275        let cuboid = ConvexPolygon::with_vertices([
276            [-1.0, -2.0].into(),
277            [1.0, -2.0].into(),
278            [1.0, 2.0].into(),
279            [-1.0, 2.0].into(),
280        ])
281        .expect("hard-coded vertices form a polygon");
282
283        assert_relative_eq!(
284            cuboid.support_mapping(&Cartesian::from([1.0, 0.1])),
285            [1.0, 2.0].into()
286        );
287        assert_relative_eq!(
288            cuboid.support_mapping(&Cartesian::from([1.0, -0.1])),
289            [1.0, -2.0].into()
290        );
291        assert_relative_eq!(
292            cuboid.support_mapping(&Cartesian::from([-0.1, 1.0])),
293            [-1.0, 2.0].into()
294        );
295        assert_relative_eq!(
296            cuboid.support_mapping(&Cartesian::from([-0.1, -1.0])),
297            [-1.0, -2.0].into()
298        );
299    }
300
301    #[test]
302    fn support_mapping_3d() {
303        let cuboid = ConvexPolyhedron::with_vertices([
304            [-1.0, -2.0, 3.0].into(),
305            [1.0, -2.0, 3.0].into(),
306            [1.0, 2.0, 3.0].into(),
307            [-1.0, 2.0, 3.0].into(),
308            [-1.0, -2.0, -3.0].into(),
309            [1.0, -2.0, -3.0].into(),
310            [1.0, 2.0, -3.0].into(),
311            [-1.0, 2.0, -3.0].into(),
312        ])
313        .expect("hard-coded vertices form a polygon");
314
315        assert_relative_eq!(
316            cuboid.support_mapping(&Cartesian::from([1.0, 0.1, 0.1])),
317            [1.0, 2.0, 3.0].into()
318        );
319        assert_relative_eq!(
320            cuboid.support_mapping(&Cartesian::from([1.0, 0.1, -0.1])),
321            [1.0, 2.0, -3.0].into()
322        );
323        assert_relative_eq!(
324            cuboid.support_mapping(&Cartesian::from([1.0, -0.1, 0.1])),
325            [1.0, -2.0, 3.0].into()
326        );
327        assert_relative_eq!(
328            cuboid.support_mapping(&Cartesian::from([1.0, -0.1, -0.1])),
329            [1.0, -2.0, -3.0].into()
330        );
331        assert_relative_eq!(
332            cuboid.support_mapping(&Cartesian::from([-1.0, 0.1, 0.1])),
333            [-1.0, 2.0, 3.0].into()
334        );
335        assert_relative_eq!(
336            cuboid.support_mapping(&Cartesian::from([-1.0, 0.1, -0.1])),
337            [-1.0, 2.0, -3.0].into()
338        );
339        assert_relative_eq!(
340            cuboid.support_mapping(&Cartesian::from([-1.0, -0.1, 0.1])),
341            [-1.0, -2.0, 3.0].into()
342        );
343        assert_relative_eq!(
344            cuboid.support_mapping(&Cartesian::from([-1.0, -0.1, -0.1])),
345            [-1.0, -2.0, -3.0].into()
346        );
347    }
348
349    // ConvexPolygon tests from hoomd-blue's test_convex_polygon.cc
350
351    #[fixture]
352    fn square() -> Convex<ConvexPolygon> {
353        Convex(
354            ConvexPolygon::with_vertices([
355                [-0.5, -0.5].into(),
356                [0.5, -0.5].into(),
357                [0.5, 0.5].into(),
358                [-0.5, 0.5].into(),
359            ])
360            .expect("hard-coded vertices form a valid polygon"),
361        )
362    }
363
364    #[fixture]
365    fn triangle() -> Convex<ConvexPolygon> {
366        Convex(
367            ConvexPolygon::with_vertices([
368                [-0.5, -0.5].into(),
369                [0.5, -0.5].into(),
370                [0.5, 0.5].into(),
371            ])
372            .expect("hard-coded vertices form a valid polygon"),
373        )
374    }
375
376    #[rstest]
377    fn square_no_rot(square: Convex<ConvexPolygon>) {
378        let a = Angle::identity();
379        assert!(!square.intersects_at(&square, &[10.0, 0.0].into(), &a));
380        assert!(!square.intersects_at(&square, &[-10.0, 0.0].into(), &a));
381
382        assert!(!square.intersects_at(&square, &[1.1, 0.0].into(), &a));
383        assert!(!square.intersects_at(&square, &[-1.1, 0.0].into(), &a));
384        assert!(!square.intersects_at(&square, &[0.0, 1.1].into(), &a));
385        assert!(!square.intersects_at(&square, &[0.0, -1.1].into(), &a));
386
387        assert!(square.intersects_at(&square, &[0.9, 0.2].into(), &a));
388        assert!(square.intersects_at(&square, &[-0.9, 0.2].into(), &a));
389        assert!(square.intersects_at(&square, &[-0.2, 0.9].into(), &a));
390        assert!(square.intersects_at(&square, &[-0.2, -0.9].into(), &a));
391
392        assert!(square.intersects_at(&square, &[1.0, 0.2].into(), &a));
393    }
394
395    #[rstest]
396    fn square_rot(square: Convex<ConvexPolygon>) {
397        let a = Angle::from(PI / 4.0);
398
399        assert!(!square.intersects_at(&square, &[10.0, 0.0].into(), &a));
400        assert!(!square.intersects_at(&square, &[-10.0, 0.0].into(), &a));
401
402        assert!(!square.intersects_at(&square, &[1.3, 0.0].into(), &a));
403        assert!(!square.intersects_at(&square, &[-1.3, 0.0].into(), &a));
404        assert!(!square.intersects_at(&square, &[0.0, 1.3].into(), &a));
405        assert!(!square.intersects_at(&square, &[0.0, -1.3].into(), &a));
406
407        assert!(!square.intersects_at(&square, &[1.3, 0.2].into(), &a));
408        assert!(!square.intersects_at(&square, &[-1.3, 0.2].into(), &a));
409        assert!(!square.intersects_at(&square, &[-0.2, 1.3].into(), &a));
410        assert!(!square.intersects_at(&square, &[-0.2, -1.3].into(), &a));
411
412        assert!(square.intersects_at(&square, &[1.2, 0.2].into(), &a));
413        assert!(square.intersects_at(&square, &[-1.2, 0.2].into(), &a));
414        assert!(square.intersects_at(&square, &[-0.2, 1.2].into(), &a));
415        assert!(square.intersects_at(&square, &[-0.2, -1.2].into(), &a));
416    }
417
418    fn test_overlap<A, B, R, const N: usize>(
419        r_ab: Cartesian<N>,
420        a: &A,
421        b: &B,
422        o_a: R,
423        o_b: &R,
424    ) -> bool
425    where
426        R: Rotation + Rotate<Cartesian<N>>,
427        A: IntersectsAt<B, Cartesian<N>, R>,
428    {
429        let r_a_inverted = o_a.inverted();
430        let v_ij = r_a_inverted.rotate(&r_ab);
431        let o_ij = o_b.combine(&r_a_inverted);
432        a.intersects_at(b, &v_ij, &o_ij)
433    }
434
435    fn assert_symmetric_overlap<A, B, R, const N: usize>(
436        r_ab: Cartesian<N>,
437        a: &A,
438        b: &B,
439        r_a: R,
440        r_b: R,
441        expected: bool,
442    ) where
443        R: Rotation + Rotate<Cartesian<N>>,
444        A: IntersectsAt<B, Cartesian<N>, R>,
445        B: IntersectsAt<A, Cartesian<N>, R>,
446    {
447        assert_eq!(test_overlap(r_ab, a, b, r_a, &r_b), expected);
448        assert_eq!(test_overlap(-r_ab, b, a, r_b, &r_a), expected);
449    }
450
451    #[rstest]
452    fn square_triangle(square: Convex<ConvexPolygon>, triangle: Convex<ConvexPolygon>) {
453        let r_square = Angle::from(-PI / 4.0);
454        let r_triangle = Angle::from(PI);
455
456        assert_symmetric_overlap(
457            [10.0, 0.0].into(),
458            &square,
459            &triangle,
460            r_square,
461            r_triangle,
462            false,
463        );
464
465        assert_symmetric_overlap(
466            [1.3, 0.0].into(),
467            &square,
468            &triangle,
469            r_square,
470            r_triangle,
471            false,
472        );
473
474        assert_symmetric_overlap(
475            [-1.3, 0.0].into(),
476            &square,
477            &triangle,
478            r_square,
479            r_triangle,
480            false,
481        );
482
483        assert_symmetric_overlap(
484            [0.0, 1.3].into(),
485            &square,
486            &triangle,
487            r_square,
488            r_triangle,
489            false,
490        );
491
492        assert_symmetric_overlap(
493            [0.0, -1.3].into(),
494            &square,
495            &triangle,
496            r_square,
497            r_triangle,
498            false,
499        );
500
501        assert_symmetric_overlap(
502            [1.2, 0.2].into(),
503            &square,
504            &triangle,
505            r_square,
506            r_triangle,
507            true,
508        );
509
510        assert_symmetric_overlap(
511            [-0.7, -0.2].into(),
512            &square,
513            &triangle,
514            r_square,
515            r_triangle,
516            true,
517        );
518
519        assert_symmetric_overlap(
520            [0.4, 1.1].into(),
521            &square,
522            &triangle,
523            r_square,
524            r_triangle,
525            true,
526        );
527
528        assert_symmetric_overlap(
529            [-0.2, -1.2].into(),
530            &square,
531            &triangle,
532            r_square,
533            r_triangle,
534            true,
535        );
536    }
537
538    #[fixture]
539    fn octahedron() -> Convex<ConvexPolyhedron> {
540        Convex(
541            ConvexPolyhedron::with_vertices([
542                [-0.5, -0.5, 0.0].into(),
543                [0.5, -0.5, 0.0].into(),
544                [0.5, 0.5, 0.0].into(),
545                [-0.5, 0.5, 0.0].into(),
546                [0.0, 0.0, FRAC_1_SQRT_2].into(),
547                [0.0, 0.0, -FRAC_1_SQRT_2].into(),
548            ])
549            .expect("hard-coded vertices form a valid polyhedron"),
550        )
551    }
552
553    #[fixture]
554    fn cube() -> Convex<ConvexPolyhedron> {
555        Convex(
556            ConvexPolyhedron::with_vertices([
557                [-0.5, -0.5, -0.5].into(),
558                [0.5, -0.5, -0.5].into(),
559                [0.5, 0.5, -0.5].into(),
560                [-0.5, 0.5, -0.5].into(),
561                [-0.5, -0.5, 0.5].into(),
562                [0.5, -0.5, 0.5].into(),
563                [0.5, 0.5, 0.5].into(),
564                [-0.5, 0.5, 0.5].into(),
565            ])
566            .expect("hard-coded vertices form a valid polyhedron"),
567        )
568    }
569
570    #[rstest]
571    fn overlap_octahedron_no_rot(octahedron: Convex<ConvexPolyhedron>) {
572        let q = Versor::identity();
573
574        assert_symmetric_overlap([0.0, 0.0, 0.0].into(), &octahedron, &octahedron, q, q, true);
575
576        assert_symmetric_overlap(
577            [10.0, 0.0, 0.0].into(),
578            &octahedron,
579            &octahedron,
580            q,
581            q,
582            false,
583        );
584
585        assert_symmetric_overlap(
586            [1.1, 0.0, 0.0].into(),
587            &octahedron,
588            &octahedron,
589            q,
590            q,
591            false,
592        );
593
594        assert_symmetric_overlap(
595            [0.0, 1.1, 0.0].into(),
596            &octahedron,
597            &octahedron,
598            q,
599            q,
600            false,
601        );
602
603        assert_symmetric_overlap(
604            [1.1, 0.2, 0.0].into(),
605            &octahedron,
606            &octahedron,
607            q,
608            q,
609            false,
610        );
611
612        assert_symmetric_overlap(
613            [-1.1, 0.2, 0.0].into(),
614            &octahedron,
615            &octahedron,
616            q,
617            q,
618            false,
619        );
620
621        assert_symmetric_overlap(
622            [-0.2, 1.1, 0.0].into(),
623            &octahedron,
624            &octahedron,
625            q,
626            q,
627            false,
628        );
629
630        assert_symmetric_overlap(
631            [-0.2, -1.1, 0.0].into(),
632            &octahedron,
633            &octahedron,
634            q,
635            q,
636            false,
637        );
638
639        assert_symmetric_overlap([0.9, 0.2, 0.0].into(), &octahedron, &octahedron, q, q, true);
640
641        assert_symmetric_overlap(
642            [-0.9, 0.2, 0.0].into(),
643            &octahedron,
644            &octahedron,
645            q,
646            q,
647            true,
648        );
649
650        assert_symmetric_overlap(
651            [-0.2, 0.9, 0.0].into(),
652            &octahedron,
653            &octahedron,
654            q,
655            q,
656            true,
657        );
658
659        assert_symmetric_overlap(
660            [-0.2, -0.9, 0.0].into(),
661            &octahedron,
662            &octahedron,
663            q,
664            q,
665            true,
666        );
667
668        assert_symmetric_overlap([1.0, 0.2, 0.0].into(), &octahedron, &octahedron, q, q, true);
669    }
670
671    #[rstest]
672    fn overlap_cube_no_rot(cube: Convex<ConvexPolyhedron>) {
673        let q = Versor::identity();
674
675        assert_symmetric_overlap([0.0, 0.0, 0.0].into(), &cube, &cube, q, q, true);
676        assert_symmetric_overlap([10.0, 0.0, 0.0].into(), &cube, &cube, q, q, false);
677
678        assert_symmetric_overlap([1.1, 0.0, 0.0].into(), &cube, &cube, q, q, false);
679        assert_symmetric_overlap([0.0, 1.1, 0.0].into(), &cube, &cube, q, q, false);
680        assert_symmetric_overlap([0.0, 0.0, 1.1].into(), &cube, &cube, q, q, false);
681        assert_symmetric_overlap([1.1, 0.2, 0.0].into(), &cube, &cube, q, q, false);
682
683        assert_symmetric_overlap([-1.1, 0.2, 0.0].into(), &cube, &cube, q, q, false);
684        assert_symmetric_overlap([-0.2, 1.1, 0.0].into(), &cube, &cube, q, q, false);
685        assert_symmetric_overlap([-0.2, -1.1, 0.0].into(), &cube, &cube, q, q, false);
686
687        assert_symmetric_overlap([0.9, 0.2, 0.0].into(), &cube, &cube, q, q, true);
688        assert_symmetric_overlap([-0.9, 0.2, 0.0].into(), &cube, &cube, q, q, true);
689        assert_symmetric_overlap([-0.2, 0.9, 0.0].into(), &cube, &cube, q, q, true);
690        assert_symmetric_overlap([-0.2, -0.9, 0.0].into(), &cube, &cube, q, q, true);
691
692        assert_symmetric_overlap([0.2, 0.0, 0.0].into(), &cube, &cube, q, q, true);
693        assert_symmetric_overlap([0.2, 0.00001, 0.00001].into(), &cube, &cube, q, q, true);
694        assert_symmetric_overlap([0.1, 0.2, 0.1].into(), &cube, &cube, q, q, true);
695        assert_symmetric_overlap([1.0, 0.2, 0.0].into(), &cube, &cube, q, q, true);
696    }
697
698    #[rstest]
699    fn overlap_cube_rot1(cube: Convex<ConvexPolyhedron>) {
700        let q_a = Versor::identity();
701        let q_b = Versor::from_axis_angle(
702            [0.0, 0.0, 1.0]
703                .try_into()
704                .expect("hard-coded vector is non-zero"),
705            PI / 4.0,
706        );
707
708        assert_symmetric_overlap([10.0, 0.0, 0.0].into(), &cube, &cube, q_a, q_b, false);
709
710        assert_symmetric_overlap([1.3, 0.0, 0.0].into(), &cube, &cube, q_a, q_b, false);
711        assert_symmetric_overlap([-1.3, 0.0, 0.0].into(), &cube, &cube, q_a, q_b, false);
712        assert_symmetric_overlap([0.0, 1.3, 0.0].into(), &cube, &cube, q_a, q_b, false);
713        assert_symmetric_overlap([0.0, -1.3, 0.0].into(), &cube, &cube, q_a, q_b, false);
714
715        assert_symmetric_overlap([1.3, 0.2, 0.0].into(), &cube, &cube, q_a, q_b, false);
716        assert_symmetric_overlap([-1.3, 0.2, 0.0].into(), &cube, &cube, q_a, q_b, false);
717        assert_symmetric_overlap([-0.2, 1.3, 0.0].into(), &cube, &cube, q_a, q_b, false);
718        assert_symmetric_overlap([-0.2, -1.3, 0.0].into(), &cube, &cube, q_a, q_b, false);
719
720        assert_symmetric_overlap([1.2, 0.2, 0.0].into(), &cube, &cube, q_a, q_b, true);
721        assert_symmetric_overlap([-1.2, 0.2, 0.0].into(), &cube, &cube, q_a, q_b, true);
722        assert_symmetric_overlap([-0.2, 1.2, 0.0].into(), &cube, &cube, q_a, q_b, true);
723        assert_symmetric_overlap([-0.2, -1.2, 0.0].into(), &cube, &cube, q_a, q_b, true);
724    }
725
726    #[rstest]
727    fn overlap_cube_rot3(cube: Convex<ConvexPolyhedron>) {
728        let q_a = Versor::identity();
729        let q1 = Versor::from_axis_angle(
730            [1.0, 0.0, 0.0]
731                .try_into()
732                .expect("hard-coded vector is non-zero"),
733            PI / 4.0,
734        );
735        let q2 = Versor::from_axis_angle(
736            [0.0, 0.0, 1.0]
737                .try_into()
738                .expect("hard-coded vector is non-zero"),
739            PI / 4.0,
740        );
741        let q_b = q2.combine(&q1);
742
743        assert_symmetric_overlap([10.0, 0.0, 0.0].into(), &cube, &cube, q_a, q_b, false);
744
745        assert_symmetric_overlap([1.4, 0.0, 0.0].into(), &cube, &cube, q_a, q_b, false);
746        assert_symmetric_overlap([-1.4, 0.0, 0.0].into(), &cube, &cube, q_a, q_b, false);
747        assert_symmetric_overlap([0.0, 1.4, 0.0].into(), &cube, &cube, q_a, q_b, false);
748        assert_symmetric_overlap([0.0, -1.4, 0.0].into(), &cube, &cube, q_a, q_b, false);
749
750        assert_symmetric_overlap([1.4, 0.2, 0.0].into(), &cube, &cube, q_a, q_b, false);
751        assert_symmetric_overlap([-1.4, 0.2, 0.0].into(), &cube, &cube, q_a, q_b, false);
752        assert_symmetric_overlap([-0.2, 1.4, 0.0].into(), &cube, &cube, q_a, q_b, false);
753        assert_symmetric_overlap([-0.2, -1.4, 0.0].into(), &cube, &cube, q_a, q_b, false);
754
755        assert_symmetric_overlap([0.0, 1.2, 0.0].into(), &cube, &cube, q_a, q_b, true);
756        assert_symmetric_overlap([0.0, 1.2, 0.1].into(), &cube, &cube, q_a, q_b, true);
757        assert_symmetric_overlap([0.1, 1.2, 0.1].into(), &cube, &cube, q_a, q_b, true);
758        assert_symmetric_overlap([1.2, 0.0, 0.0].into(), &cube, &cube, q_a, q_b, true);
759        assert_symmetric_overlap([1.2, 0.1, 0.0].into(), &cube, &cube, q_a, q_b, true);
760        assert_symmetric_overlap([1.2, 0.1, 0.1].into(), &cube, &cube, q_a, q_b, true);
761
762        assert_symmetric_overlap([-0.9, 0.9, 0.0].into(), &cube, &cube, q_a, q_b, true);
763        assert_symmetric_overlap([-0.9, 0.899, 0.001].into(), &cube, &cube, q_a, q_b, true);
764        assert_symmetric_overlap([0.9, -0.9, 0.0].into(), &cube, &cube, q_a, q_b, true);
765        assert_symmetric_overlap([-0.9, 0.9, 0.1].into(), &cube, &cube, q_a, q_b, true);
766    }
767}