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