hoomd_geometry/shape/
convex_surface_mesh_2d.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//! Convex polygon represented by vertices and edges.
5
6use std::cmp::Ordering;
7
8use itertools::Itertools;
9use serde::{Deserialize, Serialize};
10
11use crate::{
12    BoundingSphereRadius, Error, IntersectsAt, IntersectsAtGlobal, IsPointInside, SupportMapping,
13    Volume, shape::ConvexPolytope,
14};
15use hoomd_utility::valid::PositiveReal;
16use hoomd_vector::{Cartesian, InnerProduct, Metric, Rotate, Rotation, RotationMatrix};
17
18/// The vertices and edges that make up a convex polygon.
19///
20/// [`ConvexPolytope::<2>`] and [`ConvexSurfaceMesh2d`] can both represent
21/// 2d convex polygons. The first is defined *implicitly* as the convex hull
22/// of a set of points. It stores the given point set without any modification,
23/// and can therefore be constructed quickly. The *implicit* convex hull is
24/// formed by [`SupportMapping`] during intersection tests of
25/// `Convex(ConvexPolygon)` with other `Convex(_)` types.
26///
27/// In contrast, [`ConvexSurfaceMesh2d`] *explicitly* computes the convex hull
28/// on construction. After construction, the [`vertices`] of the shape include
29/// only the points on the convex hull in a counter-clockwise order.
30/// Using this representation, [`ConvexSurfaceMesh2d`] is able to provide
31/// implementations of [`Volume`], [`IsPointInside`], and [`IntersectsAt`].
32/// The native `ConvexSurfaceMesh2d`--`ConvexSurfaceMesh2d` intersection test
33/// is much faster than the generic `Convex(ConvexPolygon)` intersection test.
34///
35/// [`vertices`]: Self::vertices
36///
37/// # Examples
38///
39/// Construction:
40/// ```
41/// use hoomd_geometry::shape::ConvexSurfaceMesh2d;
42///
43/// # fn main() -> Result<(), hoomd_geometry::Error> {
44/// let triangle = ConvexSurfaceMesh2d::from_point_set([
45///     [1.0, -1.0].into(),
46///     [-1.0, -1.0].into(),
47///     [0.0, 1.0].into(),
48/// ])?;
49/// # Ok(())
50/// # }
51/// ```
52///
53/// Intersection tests:
54/// ```
55/// use hoomd_geometry::{Convex, IntersectsAt, shape::ConvexSurfaceMesh2d};
56/// use hoomd_vector::{Angle, Cartesian};
57/// use std::f64::consts::PI;
58///
59/// # fn main() -> Result<(), hoomd_geometry::Error> {
60/// let rectangle = ConvexSurfaceMesh2d::from_point_set([
61///     [-2.0, -1.0].into(),
62///     [2.0, -1.0].into(),
63///     [2.0, 1.0].into(),
64///     [-2.0, 1.0].into(),
65/// ])?;
66///
67/// assert!(!rectangle.intersects_at(
68///     &rectangle,
69///     &[0.0, 2.1].into(),
70///     &Angle::default()
71/// ));
72/// assert!(rectangle.intersects_at(
73///     &rectangle,
74///     &[0.0, 2.1].into(),
75///     &Angle::from(PI / 2.0)
76/// ));
77/// # Ok(())
78/// # }
79/// ```
80#[derive(Clone, Debug, PartialEq, Serialize, Deserialize)]
81pub struct ConvexSurfaceMesh2d {
82    /// The vertices of the polygon in counter-clockwise order.
83    vertices: Vec<Cartesian<2>>,
84    /// The radius of a bounding sphere of the geometry.
85    bounding_radius: PositiveReal,
86}
87
88/// Find the lowest, leftmost point from a slice of Cartesian vectors.
89#[inline]
90fn find_lowest_leftmost(vertices: &[Cartesian<2>]) -> Option<usize> {
91    vertices.iter().position_min_by(|a, b| {
92        // Compare y-coordinates, then x.
93        a[1].total_cmp(&b[1]).then(a[0].total_cmp(&b[0]))
94    })
95}
96
97/// Get the key for a lexicographic order of points with respect to an anchor.
98#[inline]
99fn get_graham_key(p: Cartesian<2>, anchor: Cartesian<2>) -> (f64, f64) {
100    let diff = p - anchor;
101    (f64::atan2(diff[1], diff[0]), diff.dot(&diff))
102}
103
104/// Determines whether a point `test` is to the left, right, or collinear with `edge`.
105///
106/// # Warning
107///
108/// This predicate is not robust: points very close to the line may be misclassified
109/// due to floating-point precision limits. For all practical inputs, this will not
110/// result in issues.
111///
112/// # Note
113///
114/// This formulation (often referred to as the shoelace formula) guarantees
115/// **cyclic invariance**, or the property that the orientation sign is identical
116/// for any ordering of the three points `(e0, e1, t)`, `(e1, t, e0)`, and
117/// `(t, e0, e1)`. As a result, the result is antisymmetric about the edge `e`
118/// such that `p == -p'` for any `p'` reflected over `e`.
119///
120/// These properties do *not* prevent misclassification of points near the line—they
121/// only ensure consistent behavior for related inputs.
122///
123/// **Source:** [Robust Arithmetic in Computational Geometry](https://observablehq.com/@mourner/non-robust-arithmetic-as-art)
124#[inline]
125fn predicate_orient2d((p, q): (Cartesian<2>, Cartesian<2>), test: Cartesian<2>) -> i64 {
126    let orientation = (p[0] * q[1] - p[1] * q[0])
127        + (q[0] * test[1] - q[1] * test[0])
128        + (test[0] * p[1] - test[1] * p[0]);
129
130    match orientation.total_cmp(&0.0) {
131        Ordering::Greater => 1,
132        Ordering::Less => -1,
133        Ordering::Equal => 0,
134    }
135}
136
137impl ConvexSurfaceMesh2d {
138    /// Create an convex surface mesh from the convex hull of the given set of points.
139    ///
140    /// The resulting shape contains a subset of the given points, including
141    /// only the non-degenerate points on the convex hull. The result's vertices
142    /// are arranged in a counter-clockwise order.
143    ///
144    /// # Errors
145    ///
146    /// * [`Error::DegeneratePolytope`] when there are fewer than 3 non-collinear points.
147    ///
148    /// # Example
149    /// ```
150    /// use hoomd_geometry::shape::ConvexSurfaceMesh2d;
151    ///
152    /// # fn main() -> Result<(), hoomd_geometry::Error> {
153    /// let equilateral_triangle = ConvexSurfaceMesh2d::from_point_set([
154    ///     [1.0, 0.0].into(),
155    ///     [0.5, f64::sqrt(3.0) / 2.0].into(),
156    ///     [-0.5, f64::sqrt(3.0) / 2.0].into(),
157    /// ])?;
158    /// # Ok(())
159    /// # }
160    /// ```
161    #[inline]
162    pub fn from_point_set<I>(points: I) -> Result<Self, Error>
163    where
164        I: IntoIterator<Item = Cartesian<2>>,
165    {
166        let vertices = Self::construct_convex_hull(points.into_iter().collect())?;
167
168        Ok(Self {
169            bounding_radius: ConvexPolytope::<2>::bounding_radius(&vertices),
170            vertices,
171        })
172    }
173
174    /// Compute the convex hull of points in two dimensions.
175    ///
176    /// The resulting vector contains a subset of the given points, including
177    /// only the non-degenerate points on the convex hull. The output vertices
178    /// are arranged in a counter-clockwise order.
179    ///
180    /// # Errors
181    ///
182    /// Returns [`Error`] if the input vertices do not form a convex body with 3 or more points.
183    ///
184    /// [`Error`]: enum@Error
185    ///
186    /// # Example
187    /// ```
188    /// use hoomd_geometry::shape::ConvexSurfaceMesh2d;
189    ///
190    /// # fn main() -> Result<(), hoomd_geometry::Error> {
191    /// let points = vec![
192    ///     [1.0, 0.0].into(),
193    ///     [0.5, f64::sqrt(3.0) / 2.0].into(),
194    ///     [-0.5, f64::sqrt(3.0) / 2.0].into(),
195    /// ];
196    ///
197    /// let hull_vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)?;
198    /// # Ok(())
199    /// # }
200    /// ```
201    #[inline]
202    pub fn construct_convex_hull(
203        mut points: Vec<Cartesian<2>>,
204    ) -> Result<Vec<Cartesian<2>>, Error> {
205        // No need to try and triangulate if the hull is degenerate
206        if points.len() < 3 {
207            return Err(Error::DegeneratePolytope);
208        }
209
210        let anchor_idx = find_lowest_leftmost(&points).ok_or(Error::DegeneratePolytope)?;
211
212        // Move the anchor to the front of the list of vertices, as it is always in the hull
213        points.swap(0, anchor_idx);
214        let anchor = points[0];
215
216        // Sort the remainder of the slice in-place
217        points[1..].sort_unstable_by(|&a, &b| {
218            let (a0, a1) = get_graham_key(a, anchor);
219            let (b0, b1) = get_graham_key(b, anchor);
220            a0.total_cmp(&b0).then(a1.total_cmp(&b1))
221        });
222
223        // Now vertices[..2] is an edge on the hull. Initialize counters for the hull length
224        // and number of vertices on the hull
225        let mut n_vertices_on_hull = 2;
226        let mut next_candidate = 2;
227
228        // Repeat until all interior points are gone
229        while next_candidate < points.len() {
230            let c = points[next_candidate];
231            while n_vertices_on_hull >= 2 {
232                let p = points[n_vertices_on_hull - 2];
233                let n = points[n_vertices_on_hull - 1];
234
235                if predicate_orient2d((p, n), c) <= 0 {
236                    // Point n is inside the hull, remove it by shrinking the hull
237                    n_vertices_on_hull -= 1;
238                } else {
239                    break;
240                }
241            }
242            // Swap the vertex c onto the end of the hull, extending it by one
243            points.swap(next_candidate, n_vertices_on_hull);
244            n_vertices_on_hull += 1;
245            next_candidate += 1;
246        }
247
248        points.truncate(n_vertices_on_hull);
249        if points.len() >= 3 {
250            Ok(points)
251        } else {
252            Err(Error::DegeneratePolytope)
253        }
254    }
255
256    /// The vertices of the convex polygon in counter-clockwise order.
257    #[inline]
258    #[must_use]
259    pub fn vertices(&self) -> &[Cartesian<2>] {
260        &self.vertices
261    }
262}
263
264impl SupportMapping<Cartesian<2>> for ConvexSurfaceMesh2d {
265    /// Find the point on a shape that is the furthest in a given direction.
266    ///
267    /// [`ConvexSurfaceMesh2d`] implements [`SupportMapping`] to enable
268    /// intersection tests between mixed convex types.
269    ///
270    /// # Example
271    ///
272    /// ```
273    /// use hoomd_geometry::{
274    ///     Convex, IntersectsAt,
275    ///     shape::{Circle, ConvexSurfaceMesh2d, Sphero},
276    /// };
277    /// use hoomd_vector::{Angle, Cartesian};
278    /// use std::f64::consts::PI;
279    ///
280    /// # fn main() -> Result<(), Box<dyn std::error::Error>> {
281    /// let circle = Convex(Circle {
282    ///     radius: 0.5.try_into()?,
283    /// });
284    /// let rectangle = ConvexSurfaceMesh2d::from_point_set([
285    ///     [-1.5, -1.0].into(),
286    ///     [1.5, -1.0].into(),
287    ///     [-1.5, 1.0].into(),
288    ///     [1.5, 1.0].into(),
289    /// ])?;
290    /// let rounded_rectangle = Convex(Sphero {
291    ///     shape: rectangle,
292    ///     rounding_radius: 0.5.try_into()?,
293    /// });
294    ///
295    /// assert!(rounded_rectangle.intersects_at(
296    ///     &circle,
297    ///     &[2.4, 0.0].into(),
298    ///     &Angle::default()
299    /// ));
300    /// assert!(!rounded_rectangle.intersects_at(
301    ///     &circle,
302    ///     &[0.0, 2.4].into(),
303    ///     &Angle::default()
304    /// ));
305    /// assert!(circle.intersects_at(
306    ///     &rounded_rectangle,
307    ///     &[0.0, 2.4].into(),
308    ///     &Angle::from(PI / 2.0)
309    /// ));
310    /// # Ok(())
311    /// # }
312    /// ```
313    #[inline]
314    fn support_mapping(&self, n: &Cartesian<2>) -> Cartesian<2> {
315        *self
316            .vertices
317            .iter()
318            .max_by(|a, b| {
319                a.dot(n)
320                    .partial_cmp(&b.dot(n))
321                    .unwrap_or(std::cmp::Ordering::Equal)
322            })
323            .expect("there should be at least 3 vertices")
324    }
325}
326
327impl BoundingSphereRadius for ConvexSurfaceMesh2d {
328    /// Radius of a circle that bounds the shape.
329    ///
330    /// The circle has the same local origin as the shape `self`.
331    ///
332    /// # Example
333    ///
334    /// ```
335    /// use approxim::assert_relative_eq;
336    /// use hoomd_geometry::{BoundingSphereRadius, shape::ConvexSurfaceMesh2d};
337    ///
338    /// # fn main() -> Result<(), hoomd_geometry::Error> {
339    /// let triangle = ConvexSurfaceMesh2d::from_point_set([
340    ///     [1.0, -1.0].into(),
341    ///     [-1.0, -1.0].into(),
342    ///     [0.0, 1.0].into(),
343    /// ])?;
344    ///
345    /// let bounding_radius = triangle.bounding_sphere_radius();
346    ///
347    /// assert_relative_eq!(bounding_radius.get(), 2.0_f64.sqrt());
348    /// # Ok(())
349    /// # }
350    /// ```
351    #[inline]
352    fn bounding_sphere_radius(&self) -> PositiveReal {
353        self.bounding_radius
354    }
355}
356
357impl Volume for ConvexSurfaceMesh2d {
358    /// Compute the area of the convex polygon.
359    ///
360    /// # Example
361    ///
362    /// ```
363    /// use approxim::assert_relative_eq;
364    /// use hoomd_geometry::{Volume, shape::ConvexSurfaceMesh2d};
365    ///
366    /// # fn main() -> Result<(), hoomd_geometry::Error> {
367    /// let triangle = ConvexSurfaceMesh2d::from_point_set([
368    ///     [-1.0, -1.0].into(),
369    ///     [1.0, -1.0].into(),
370    ///     [1.0, 1.0].into(),
371    /// ])?;
372    ///
373    /// let area = triangle.volume();
374    ///
375    /// assert_relative_eq!(area, 2.0);
376    /// # Ok(())
377    /// # }
378    /// ```
379    #[inline]
380    fn volume(&self) -> f64 {
381        // Compute the polygon area with the shoelace formula:
382        // https://mathworld.wolfram.com/PolygonArea.html
383        0.5 * self
384            .vertices
385            .iter()
386            .circular_tuple_windows()
387            .fold(0.0, |total, (a, b)| total + a[0] * b[1] - b[0] * a[1])
388    }
389}
390
391impl IsPointInside<Cartesian<2>> for ConvexSurfaceMesh2d {
392    /// Check if a point is inside the convex polygon.
393    ///
394    /// # Example
395    ///
396    /// ```
397    /// use hoomd_geometry::{IsPointInside, shape::ConvexSurfaceMesh2d};
398    ///
399    /// # fn main() -> Result<(), hoomd_geometry::Error> {
400    /// let triangle = ConvexSurfaceMesh2d::from_point_set([
401    ///     [1.0, -1.0].into(),
402    ///     [-1.0, -1.0].into(),
403    ///     [0.0, 1.0].into(),
404    /// ])?;
405    ///
406    /// assert!(triangle.is_point_inside(&[0.0, 0.0].into()));
407    /// assert!(triangle.is_point_inside(&[-0.9, -0.9].into()));
408    /// assert!(!triangle.is_point_inside(&[-1.5, 2.0].into()));
409    /// # Ok(())
410    /// # }
411    /// ```
412    #[inline]
413    fn is_point_inside(&self, point: &Cartesian<2>) -> bool {
414        for (a, b) in self.vertices.iter().circular_tuple_windows() {
415            let edge = *b - *a;
416            let n = -edge.perpendicular();
417
418            let v = *point - *a;
419            if v.dot(&n) > 0.0 {
420                return false;
421            }
422        }
423
424        true
425    }
426}
427
428impl<R> IntersectsAtGlobal<Self, Cartesian<2>, R> for ConvexSurfaceMesh2d
429where
430    R: Rotation + Rotate<Cartesian<2>>,
431    RotationMatrix<2>: From<R>,
432{
433    #[inline]
434    fn intersects_at_global(
435        &self,
436        other: &Self,
437        r_self: &Cartesian<2>,
438        o_self: &R,
439        r_other: &Cartesian<2>,
440        o_other: &R,
441    ) -> bool {
442        let max_separation =
443            self.bounding_sphere_radius().get() + other.bounding_sphere_radius().get();
444        if r_self.distance_squared(r_other) >= max_separation.powi(2) {
445            return false;
446        }
447
448        let (v_ij, o_ij) = hoomd_vector::pair_system_to_local(r_self, o_self, r_other, o_other);
449
450        self.intersects_at(other, &v_ij, &o_ij)
451    }
452}
453
454impl<R> IntersectsAt<Self, Cartesian<2>, R> for ConvexSurfaceMesh2d
455where
456    RotationMatrix<2>: From<R>,
457    R: Copy,
458{
459    /// Test convex polygon intersections using the separating planes method.
460    ///
461    /// When the number of vertices is small, separating planes is significantly
462    /// faster than the Xenocollide algorithm implemented for the `Convex<ConvexPolygon>`
463    /// type, even though separating planes is $` O(n^2) `$.
464    ///
465    /// # Example
466    /// ```
467    /// use hoomd_geometry::{IntersectsAt, shape::ConvexSurfaceMesh2d};
468    /// use hoomd_vector::{Angle, Cartesian};
469    /// use std::f64::consts::PI;
470    ///
471    /// # fn main() -> Result<(), hoomd_geometry::Error> {
472    /// let rectangle = ConvexSurfaceMesh2d::from_point_set([
473    ///     [-2.0, -1.0].into(),
474    ///     [2.0, -1.0].into(),
475    ///     [2.0, 1.0].into(),
476    ///     [-2.0, 1.0].into(),
477    /// ])?;
478    ///
479    /// assert!(!rectangle.intersects_at(
480    ///     &rectangle,
481    ///     &[0.0, 2.1].into(),
482    ///     &Angle::default()
483    /// ));
484    /// assert!(rectangle.intersects_at(
485    ///     &rectangle,
486    ///     &[0.0, 2.1].into(),
487    ///     &Angle::from(PI / 2.0)
488    /// ));
489    /// # Ok(())
490    /// # }
491    /// ```
492    #[inline]
493    fn intersects_at(&self, other: &Self, v_ij: &Cartesian<2>, o_ij: &R) -> bool {
494        assert!(
495            self.vertices.len() > 2 && other.vertices.len() > 2,
496            "A convex polygon must have at least 3 vertices."
497        );
498
499        let o_j = RotationMatrix::from(*o_ij);
500        if b_edge_separates(self, other, v_ij, &o_j) {
501            return false;
502        }
503
504        let o_j_inverted = o_j.inverted();
505        let v_ji = o_j_inverted.rotate(&-*v_ij);
506        if b_edge_separates(other, self, &v_ji, &o_j_inverted) {
507            return false;
508        }
509
510        true
511    }
512}
513
514/// Determine if any edge of `b` separates the points in `a` and `b`.
515#[inline]
516fn b_edge_separates(
517    a: &ConvexSurfaceMesh2d,
518    b: &ConvexSurfaceMesh2d,
519    v_ab: &Cartesian<2>,
520    o_b: &RotationMatrix<2>,
521) -> bool {
522    let mut previous = b.vertices.len() - 1;
523    for current in 0..b.vertices.len() {
524        let p = b.vertices[current];
525        let edge = p - b.vertices[previous];
526
527        let n = -edge.perpendicular();
528
529        let p_in_frame_a = o_b.rotate(&p) + *v_ab;
530        let n_in_frame_a = o_b.rotate(&n);
531
532        if is_separating(a, &p_in_frame_a, &n_in_frame_a) {
533            return true;
534        }
535
536        previous = current;
537    }
538
539    false
540}
541
542/// Determine if all of a's vertices are outside the given plane.
543#[inline]
544fn is_separating(a: &ConvexSurfaceMesh2d, p: &Cartesian<2>, n: &Cartesian<2>) -> bool {
545    // check if n dot (v[i]-p) < 0 for every vertex in the polygon
546    // distribute: (n dot v[i] - n dot p) < 0
547    let n_dot_p = n.dot(p);
548
549    for v in &a.vertices {
550        if n.dot(v) - n_dot_p <= 0.0 {
551            return false;
552        }
553    }
554
555    true
556}
557
558impl<const MAX_VERTICES: usize> TryFrom<ConvexPolytope<2, MAX_VERTICES>> for ConvexSurfaceMesh2d {
559    type Error = Error;
560
561    /// Construct the convex hull of a [`ConvexPolytope<2>`].
562    ///
563    /// # Errors
564    ///
565    /// * [`Error::DegeneratePolytope`] when there are fewer than 3 non-collinear points.
566    ///
567    /// # Example
568    ///
569    /// Invalid conversion
570    /// ```
571    /// use hoomd_geometry::shape::{ConvexPolygon, ConvexSurfaceMesh2d};
572    ///
573    /// # fn main() -> Result<(), hoomd_geometry::Error> {
574    /// let equilateral_triangle = ConvexPolygon::regular(3);
575    /// let mesh = ConvexSurfaceMesh2d::try_from(equilateral_triangle)?;
576    /// # Ok(())
577    /// # }
578    /// ```
579    #[inline]
580    fn try_from(v: ConvexPolytope<2, MAX_VERTICES>) -> Result<ConvexSurfaceMesh2d, Error> {
581        Self::from_point_set(v.vertices().iter().copied())
582    }
583}
584
585#[cfg(test)]
586mod tests {
587    use std::f64::consts::PI;
588
589    use super::*;
590    use approxim::assert_relative_eq;
591    use assert2::check;
592    use hoomd_vector::Angle;
593    use rand::{RngExt, SeedableRng, rngs::StdRng};
594    use rstest::*;
595
596    #[rstest]
597    #[case::single_point(vec![[1.0, 2.0]], 0)]
598    #[case::two_points_first_lower(vec![[1.0, 1.0], [2.0, 3.0]], 0)]
599    #[case::two_points_second_lower(vec![[1.0, 3.0], [2.0, 1.0]], 1)]
600    #[case::same_y_leftmost_wins(vec![[3.0, 1.0], [1.0, 1.0], [2.0, 1.0]], 1)]
601    #[case::same_y_negative(vec![[0.0, -5.0], [-3.0, -5.0], [2.0, -5.0]], 1)]
602    #[case::multiple_points(vec![[3.0, 5.0], [1.0, 1.0], [4.0, 2.0], [2.0, 3.0]], 1)]
603    #[case::negative_coords(vec![[0.0, 0.0], [-1.0, -1.0], [1.0, -1.0]], 1)]
604    #[case::same_y_all_negative_x(vec![[5.0, 0.0], [-10.0, 0.0], [-5.0, 0.0]], 1)]
605    #[case::lowest_is_only_point(vec![[100.0, -100.0]], 0)]
606    #[case::diagonal_tiebreak(vec![[5.0, 5.0], [4.0, 4.0], [3.0, 3.0], [2.0, 2.0], [1.0, 1.0]], 4)]
607    fn test_find_lowest_leftmost(#[case] vertices: Vec<[f64; 2]>, #[case] expected_idx: usize) {
608        let vertices: Vec<Cartesian<2>> = vertices.into_iter().map(Cartesian::from).collect();
609        let idx = find_lowest_leftmost(&vertices).expect("returned None for non-empty input");
610        assert_eq!(idx, expected_idx);
611    }
612
613    #[rstest]
614    fn test_find_lowest_leftmost_empty() {
615        let vertices: Vec<Cartesian<2>> = vec![];
616        assert_eq!(find_lowest_leftmost(&vertices), None);
617    }
618
619    #[rstest]
620    fn test_single_point_various_coords(
621        #[values([0.0, 0.0], [-1.0, -1.0], [100.0, -50.0], [f64::MIN_POSITIVE, f64::MIN_POSITIVE])]
622        coords: [f64; 2],
623    ) {
624        let vertices = vec![Cartesian::from(coords)];
625        let idx = find_lowest_leftmost(&vertices).expect("returned None for non-empty input");
626        assert_eq!(idx, 0);
627    }
628
629    #[rstest]
630    fn test_square_corners() {
631        let points: Vec<Cartesian<2>> = vec![[1.0, 1.0], [1.0, 0.0], [0.0, 0.0], [0.0, 1.0]]
632            .into_iter()
633            .map(Cartesian::from)
634            .collect();
635        let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
636            .expect("hard-coded points should lie on a convex hull");
637        assert_eq!(vertices.len(), 4);
638
639        let hull = [
640            [0.0, 0.0].into(),
641            [1.0, 0.0].into(),
642            [1.0, 1.0].into(),
643            [0.0, 1.0].into(),
644        ];
645        itertools::assert_equal(&vertices, &hull);
646    }
647
648    #[rstest]
649    fn test_square_corners_big() {
650        let points: Vec<Cartesian<2>> = vec![
651            [101.0, 101.0],
652            [101.0, 100.0],
653            [100.0, 101.0],
654            [100.0, 100.0],
655        ]
656        .into_iter()
657        .map(Cartesian::from)
658        .collect();
659        let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
660            .expect("hard-coded points should lie on a convex hull");
661        assert_eq!(vertices.len(), 4);
662
663        let hull = [
664            [100.0, 100.0].into(),
665            [101.0, 100.0].into(),
666            [101.0, 101.0].into(),
667            [100.0, 101.0].into(),
668        ];
669        itertools::assert_equal(&vertices, &hull);
670    }
671
672    #[rstest]
673    fn test_square_with_edge_points() {
674        let points: Vec<Cartesian<2>> = vec![
675            [0.0, 0.0],
676            [0.5, 0.0],
677            [1.0, 0.0], // bottom edge
678            [1.0, 0.5],
679            [1.0, 1.0], // right edge
680            [0.5, 1.0],
681            [0.0, 1.0], // top edge
682            [0.0, 0.5], // left edge
683        ]
684        .into_iter()
685        .map(Cartesian::from)
686        .collect();
687        let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
688            .expect("hard-coded points should lie on a convex hull");
689        // Hull should have 4 corners (interior edge points excluded)
690        assert_eq!(vertices.len(), 4);
691        let hull = [
692            [0.0, 0.0].into(),
693            [1.0, 0.0].into(),
694            [1.0, 1.0].into(),
695            [0.0, 1.0].into(),
696        ];
697        itertools::assert_equal(&vertices, &hull);
698    }
699
700    #[rstest]
701    fn test_square_dense_boundary() {
702        let mut pts: Vec<[f64; 2]> = Vec::new();
703        // Bottom edge
704        for i in 0..20 {
705            pts.push([f64::from(i) / 19.0, 0.0]);
706        }
707        // Right edge
708        for i in 0..20 {
709            pts.push([1.0, f64::from(i) / 19.0]);
710        }
711        // Top edge
712        for i in 0..20 {
713            pts.push([f64::from(i) / 19.0, 1.0]);
714        }
715        // Left edge
716        for i in 0..20 {
717            pts.push([0.0, f64::from(i) / 19.0]);
718        }
719        let points: Vec<Cartesian<2>> = pts.into_iter().map(Cartesian::from).collect();
720        let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
721            .expect("hard-coded points should lie on a convex hull");
722        assert_eq!(vertices.len(), 4);
723
724        let hull = [
725            [0.0, 0.0].into(),
726            [1.0, 0.0].into(),
727            [1.0, 1.0].into(),
728            [0.0, 1.0].into(),
729        ];
730        itertools::assert_equal(&vertices, &hull);
731    }
732
733    #[rstest]
734    fn test_circle_uniform() {
735        let n = 20;
736        let points: Vec<Cartesian<2>> = (0..n)
737            .map(|i| {
738                let angle = 2.0 * std::f64::consts::PI * i as f64 / n as f64;
739                Cartesian::from([angle.cos(), angle.sin()])
740            })
741            .collect();
742        let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
743            .expect("hard-coded points should lie on a convex hull");
744        // All points on circle should be in hull
745        assert_eq!(vertices.len(), n);
746    }
747
748    #[rstest]
749    fn test_circle_with_interior_points() {
750        let n_boundary = 12;
751        let mut points: Vec<Cartesian<2>> = (0..n_boundary)
752            .map(|i| {
753                let angle = 2.0 * std::f64::consts::PI * i as f64 / n_boundary as f64;
754                Cartesian::from([angle.cos(), angle.sin()])
755            })
756            .collect();
757        // Add interior points
758        points.extend([[0.0, 0.0], [0.3, 0.3], [-0.2, 0.1], [0.1, -0.4]].map(Cartesian::from));
759        let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
760            .expect("hard-coded points should lie on a convex hull");
761        // Only boundary points should be in hull
762        assert_eq!(vertices.len(), n_boundary);
763    }
764
765    #[rstest]
766    fn test_circle_partial_arc() {
767        let points: Vec<Cartesian<2>> = (0..10)
768            .map(|i| {
769                let angle = -std::f64::consts::FRAC_PI_4
770                    + (std::f64::consts::FRAC_PI_2 * f64::from(i) / 9.0);
771                Cartesian::from([angle.cos(), angle.sin()])
772            })
773            .collect();
774        let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
775            .expect("hard-coded points should lie on a convex hull");
776        // All points on partial arc should be in hull
777        assert_eq!(vertices.len(), 10);
778    }
779
780    #[rstest]
781    fn test_random_unit_square(#[values(0, 42, 64, 100_000)] seed: u64) {
782        let mut rng = StdRng::seed_from_u64(seed);
783        let original: Vec<Cartesian<2>> = (0..50)
784            .map(|_| Cartesian::from([rng.random::<f64>(), rng.random::<f64>()]))
785            .collect();
786        let points = original.clone();
787        let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
788            .expect("hard-coded points should lie on a convex hull");
789        // Hull should have at least 3 points
790        assert!(vertices.len() >= 3);
791        // All hull points should be in original set
792        for h in &vertices {
793            assert!(original.iter().any(|v| (*v - *h).norm() < 1e-10));
794        }
795    }
796
797    #[rstest]
798    fn test_random_gaussian() {
799        let mut rng = StdRng::seed_from_u64(42);
800        let points: Vec<Cartesian<2>> = (0..100)
801            .map(|_| {
802                Cartesian::from([
803                    rng.random::<f64>() * 2.0 - 1.0,
804                    rng.random::<f64>() * 2.0 - 1.0,
805                ])
806            })
807            .collect();
808        let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
809            .expect("hard-coded points should lie on a convex hull");
810        assert!(vertices.len() >= 3);
811    }
812
813    #[rstest]
814    fn test_random_deterministic_output() {
815        for _ in 0..3 {
816            let mut rng = StdRng::seed_from_u64(123);
817            let points1: Vec<Cartesian<2>> = (0..30)
818                .map(|_| Cartesian::from([rng.random::<f64>(), rng.random::<f64>()]))
819                .collect();
820            let vertices1 = ConvexSurfaceMesh2d::construct_convex_hull(points1)
821                .expect("hard-coded points should lie on a convex hull");
822
823            let mut rng = StdRng::seed_from_u64(123);
824            let points2: Vec<Cartesian<2>> = (0..30)
825                .map(|_| Cartesian::from([rng.random::<f64>(), rng.random::<f64>()]))
826                .collect();
827            let vertices2 = ConvexSurfaceMesh2d::construct_convex_hull(points2)
828                .expect("hard-coded points should lie on a convex hull");
829
830            assert_eq!(vertices1.len(), vertices2.len());
831        }
832    }
833
834    #[rstest]
835    fn test_duplicate_lowest_points() {
836        let points: Vec<Cartesian<2>> = vec![
837            [0.0, 0.0],
838            [0.0, 0.0],
839            [0.0, 0.0], // Three duplicates of lowest
840            [1.0, 0.0],
841            [1.0, 1.0],
842            [0.0, 1.0],
843        ]
844        .into_iter()
845        .map(Cartesian::from)
846        .collect();
847        let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
848            .expect("hard-coded points should lie on a convex hull");
849        assert!(vertices.len() >= 3);
850
851        let hull = [
852            [0.0, 0.0].into(),
853            [1.0, 0.0].into(),
854            [1.0, 1.0].into(),
855            [0.0, 1.0].into(),
856        ];
857        itertools::assert_equal(&vertices, &hull);
858    }
859
860    #[rstest]
861    fn test_many_duplicates_few_unique() {
862        let points: Vec<Cartesian<2>> = vec![
863            [0.0, 0.0],
864            [0.0, 0.0],
865            [0.0, 0.0],
866            [2.0, 0.0],
867            [2.0, 0.0],
868            [1.0, 2.0],
869            [1.0, 2.0],
870            [1.0, 2.0],
871        ]
872        .into_iter()
873        .map(Cartesian::from)
874        .collect();
875        let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
876            .expect("hard-coded points should lie on a convex hull");
877        // Should handle duplicates gracefully
878        assert!(vertices.len() >= 3);
879
880        let hull = [[0.0, 0.0].into(), [2.0, 0.0].into(), [1.0, 2.0].into()];
881        itertools::assert_equal(&vertices, &hull);
882    }
883
884    #[rstest]
885    fn test_collinear_bottom_edge() {
886        let points: Vec<Cartesian<2>> = vec![
887            [0.0, 0.0],
888            [0.5, 0.0],
889            [1.0, 0.0],
890            [1.5, 0.0],  // All at y=0
891            [0.75, 1.0], // Apex
892        ]
893        .into_iter()
894        .map(Cartesian::from)
895        .collect();
896        let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
897            .expect("hard-coded points should lie on a convex hull");
898        // Should pick leftmost of bottom points and include apex
899        assert!(vertices.len() >= 3);
900
901        let hull = [[0.0, 0.0].into(), [1.5, 0.0].into(), [0.75, 1.0].into()];
902        itertools::assert_equal(&vertices, &hull);
903    }
904
905    #[rstest]
906    fn test_leftmost_selected() {
907        let points: Vec<Cartesian<2>> = vec![
908            [0.5, 0.0],
909            [1.0, 0.0],
910            [0.0, 0.0], // All y=0, but [0,0] is leftmost
911            [0.5, 1.0],
912        ]
913        .into_iter()
914        .map(Cartesian::from)
915        .collect();
916        let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
917            .expect("hard-coded points should lie on a convex hull");
918        // [0, 0] should be the anchor point
919        let hull = [[0.0, 0.0].into(), [1.0, 0.0].into(), [0.5, 1.0].into()];
920        itertools::assert_equal(&vertices, &hull);
921    }
922
923    #[rstest]
924    fn test_many_bottom_points() {
925        let points: Vec<Cartesian<2>> = (0..10)
926            .map(|i| [f64::from(i) * 2.0 / 9.0, 0.0])
927            .chain([[1.0, 1.0]])
928            .map(Cartesian::from)
929            .collect();
930        let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
931            .expect("hard-coded points should lie on a convex hull");
932        // Leftmost [0,0] and rightmost [2,0] should be in hull with apex
933        assert!(vertices.len() >= 3);
934    }
935
936    #[rstest]
937    fn test_collinear_from_anchor() {
938        let points: Vec<Cartesian<2>> = vec![
939            [0.0, 0.0], // Anchor (lowest leftmost)
940            [1.0, 1.0],
941            [2.0, 2.0],
942            [3.0, 3.0], // Collinear at 45 degrees
943            [1.0, 0.0],
944            [0.0, 1.0], // Other corners
945        ]
946        .into_iter()
947        .map(Cartesian::from)
948        .collect();
949        let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
950            .expect("hard-coded points should lie on a convex hull");
951        // Should only keep furthest point in each direction
952        assert!(vertices.len() >= 3);
953
954        let hull = [
955            [0.0, 0.0].into(),
956            [1.0, 0.0].into(),
957            [3.0, 3.0].into(),
958            [0.0, 1.0].into(),
959        ];
960        itertools::assert_equal(&vertices, &hull);
961    }
962
963    #[rstest]
964    fn test_radial_collinear_multiple_directions() {
965        let points: Vec<Cartesian<2>> = vec![
966            [0.0, 0.0], // Anchor
967            [1.0, 0.0],
968            [2.0, 0.0],
969            [3.0, 0.0], // Along x-axis
970            [0.0, 1.0],
971            [0.0, 2.0],
972            [0.0, 3.0], // Along y-axis
973            [1.0, 1.0],
974            [2.0, 2.0], // Diagonal
975            [-1.0, 0.0],
976            [-2.0, 0.0], // Negative x-axis (but won't be picked due to anchor)
977        ]
978        .into_iter()
979        .map(Cartesian::from)
980        .collect();
981        let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
982            .expect("hard-coded points should lie on a convex hull");
983        // Should keep only outermost points
984        assert!(vertices.len() >= 3);
985    }
986
987    #[rstest]
988    fn test_star_pattern() {
989        let mut points: Vec<Cartesian<2>> = vec![Cartesian::from([0.0, 0.0])]; // Anchor
990        // Create points along 8 rays
991        for i in 0..8 {
992            let angle = 2.0 * std::f64::consts::PI * f64::from(i) / 8.0;
993            for r in [0.5, 1.0, 1.5] {
994                points.push(Cartesian::from([r * angle.cos(), r * angle.sin()]));
995            }
996        }
997        let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
998            .expect("hard-coded points should lie on a convex hull");
999        // Outer points should form the hull
1000        assert!(vertices.len() >= 8); // At least 8 outer points
1001    }
1002
1003    #[rstest]
1004    fn test_minimum_triangle() {
1005        let points: Vec<Cartesian<2>> = vec![[0.0, 0.0], [1.0, 0.0], [0.5, 1.0]]
1006            .into_iter()
1007            .map(Cartesian::from)
1008            .collect();
1009        let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
1010            .expect("hard-coded points should lie on a convex hull");
1011        assert_eq!(vertices.len(), 3);
1012
1013        let hull = [[0.0, 0.0].into(), [1.0, 0.0].into(), [0.5, 1.0].into()];
1014        itertools::assert_equal(&vertices, &hull);
1015    }
1016
1017    #[rstest]
1018    fn test_negative_coordinates() {
1019        let points: Vec<Cartesian<2>> = vec![[1.0, 1.0], [1.0, -1.0], [-1.0, 1.0], [-1.0, -1.0]]
1020            .into_iter()
1021            .map(Cartesian::from)
1022            .collect();
1023        let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
1024            .expect("hard-coded points should lie on a convex hull");
1025        assert_eq!(vertices.len(), 4);
1026
1027        let hull = [
1028            [-1.0, -1.0].into(),
1029            [1.0, -1.0].into(),
1030            [1.0, 1.0].into(),
1031            [-1.0, 1.0].into(),
1032        ];
1033        itertools::assert_equal(&vertices, &hull);
1034    }
1035
1036    #[rstest]
1037    fn test_large_coordinates() {
1038        let points: Vec<Cartesian<2>> = vec![[0.0, 0.0], [1e6, 0.0], [1e6, 1e6], [0.0, 1e6]]
1039            .into_iter()
1040            .map(Cartesian::from)
1041            .collect();
1042        let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
1043            .expect("hard-coded points should lie on a convex hull");
1044        assert_eq!(vertices.len(), 4);
1045
1046        let hull = [
1047            [0.0, 0.0].into(),
1048            [1e6, 0.0].into(),
1049            [1e6, 1e6].into(),
1050            [0.0, 1e6].into(),
1051        ];
1052        itertools::assert_equal(&vertices, &hull);
1053    }
1054
1055    #[rstest]
1056    fn test_degenerate() {
1057        let points: Vec<Cartesian<2>> = vec![[0.0, 0.0], [0.5, 0.5], [0.25, 0.25], [1.0, 1.0]]
1058            .into_iter()
1059            .map(Cartesian::from)
1060            .collect();
1061        let result = ConvexSurfaceMesh2d::construct_convex_hull(points);
1062        check!(result == Err(Error::DegeneratePolytope));
1063    }
1064
1065    #[test]
1066    fn support_mapping_2d() {
1067        let cuboid = ConvexSurfaceMesh2d::from_point_set([
1068            [-1.0, -2.0].into(),
1069            [1.0, -2.0].into(),
1070            [1.0, 2.0].into(),
1071            [-1.0, 2.0].into(),
1072        ])
1073        .expect("hard-coded vertices form a polygon");
1074
1075        assert_relative_eq!(
1076            cuboid.support_mapping(&Cartesian::from([1.0, 0.1])),
1077            [1.0, 2.0].into()
1078        );
1079        assert_relative_eq!(
1080            cuboid.support_mapping(&Cartesian::from([1.0, -0.1])),
1081            [1.0, -2.0].into()
1082        );
1083        assert_relative_eq!(
1084            cuboid.support_mapping(&Cartesian::from([-0.1, 1.0])),
1085            [-1.0, 2.0].into()
1086        );
1087        assert_relative_eq!(
1088            cuboid.support_mapping(&Cartesian::from([-0.1, -1.0])),
1089            [-1.0, -2.0].into()
1090        );
1091    }
1092
1093    // ConvexPolygon tests from hoomd-blue's test_convex_polygon.cc
1094
1095    #[fixture]
1096    fn square() -> ConvexSurfaceMesh2d {
1097        ConvexSurfaceMesh2d::from_point_set([
1098            [-0.5, -0.5].into(),
1099            [0.5, -0.5].into(),
1100            [0.5, 0.5].into(),
1101            [-0.5, 0.5].into(),
1102        ])
1103        .expect("hard-coded vertices form a valid polygon")
1104    }
1105
1106    #[fixture]
1107    fn triangle() -> ConvexSurfaceMesh2d {
1108        ConvexSurfaceMesh2d::from_point_set([
1109            [-0.5, -0.5].into(),
1110            [0.5, -0.5].into(),
1111            [0.5, 0.5].into(),
1112        ])
1113        .expect("hard-coded vertices form a valid polygon")
1114    }
1115
1116    #[rstest]
1117    fn square_no_rot(square: ConvexSurfaceMesh2d) {
1118        let a = Angle::identity();
1119        assert!(!square.intersects_at(&square, &[10.0, 0.0].into(), &a));
1120        assert!(!square.intersects_at(&square, &[-10.0, 0.0].into(), &a));
1121
1122        assert!(!square.intersects_at(&square, &[1.1, 0.0].into(), &a));
1123        assert!(!square.intersects_at(&square, &[-1.1, 0.0].into(), &a));
1124        assert!(!square.intersects_at(&square, &[0.0, 1.1].into(), &a));
1125        assert!(!square.intersects_at(&square, &[0.0, -1.1].into(), &a));
1126
1127        assert!(square.intersects_at(&square, &[0.9, 0.2].into(), &a));
1128        assert!(square.intersects_at(&square, &[-0.9, 0.2].into(), &a));
1129        assert!(square.intersects_at(&square, &[-0.2, 0.9].into(), &a));
1130        assert!(square.intersects_at(&square, &[-0.2, -0.9].into(), &a));
1131
1132        assert!(square.intersects_at(&square, &[1.0, 0.2].into(), &a));
1133    }
1134
1135    #[rstest]
1136    fn square_rot(square: ConvexSurfaceMesh2d) {
1137        let a = Angle::from(PI / 4.0);
1138
1139        assert!(!square.intersects_at(&square, &[10.0, 0.0].into(), &a));
1140        assert!(!square.intersects_at(&square, &[-10.0, 0.0].into(), &a));
1141
1142        assert!(!square.intersects_at(&square, &[1.3, 0.0].into(), &a));
1143        assert!(!square.intersects_at(&square, &[-1.3, 0.0].into(), &a));
1144        assert!(!square.intersects_at(&square, &[0.0, 1.3].into(), &a));
1145        assert!(!square.intersects_at(&square, &[0.0, -1.3].into(), &a));
1146
1147        assert!(!square.intersects_at(&square, &[1.3, 0.2].into(), &a));
1148        assert!(!square.intersects_at(&square, &[-1.3, 0.2].into(), &a));
1149        assert!(!square.intersects_at(&square, &[-0.2, 1.3].into(), &a));
1150        assert!(!square.intersects_at(&square, &[-0.2, -1.3].into(), &a));
1151
1152        assert!(square.intersects_at(&square, &[1.2, 0.2].into(), &a));
1153        assert!(square.intersects_at(&square, &[-1.2, 0.2].into(), &a));
1154        assert!(square.intersects_at(&square, &[-0.2, 1.2].into(), &a));
1155        assert!(square.intersects_at(&square, &[-0.2, -1.2].into(), &a));
1156    }
1157
1158    fn test_overlap<A, B, R, const N: usize>(
1159        r_ab: Cartesian<N>,
1160        a: &A,
1161        b: &B,
1162        o_a: R,
1163        o_b: &R,
1164    ) -> bool
1165    where
1166        R: Rotation + Rotate<Cartesian<N>>,
1167        A: IntersectsAt<B, Cartesian<N>, R>,
1168    {
1169        let r_a_inverted = o_a.inverted();
1170        let v_ij = r_a_inverted.rotate(&r_ab);
1171        let o_ij = o_b.combine(&r_a_inverted);
1172        a.intersects_at(b, &v_ij, &o_ij)
1173    }
1174
1175    fn assert_symmetric_overlap<A, B, R, const N: usize>(
1176        r_ab: Cartesian<N>,
1177        a: &A,
1178        b: &B,
1179        r_a: R,
1180        r_b: R,
1181        expected: bool,
1182    ) where
1183        R: Rotation + Rotate<Cartesian<N>>,
1184        A: IntersectsAt<B, Cartesian<N>, R>,
1185        B: IntersectsAt<A, Cartesian<N>, R>,
1186    {
1187        assert_eq!(test_overlap(r_ab, a, b, r_a, &r_b), expected);
1188        assert_eq!(test_overlap(-r_ab, b, a, r_b, &r_a), expected);
1189    }
1190
1191    #[rstest]
1192    fn square_triangle(square: ConvexSurfaceMesh2d, triangle: ConvexSurfaceMesh2d) {
1193        let r_square = Angle::from(-PI / 4.0);
1194        let r_triangle = Angle::from(PI);
1195
1196        assert_symmetric_overlap(
1197            [10.0, 0.0].into(),
1198            &square,
1199            &triangle,
1200            r_square,
1201            r_triangle,
1202            false,
1203        );
1204
1205        assert_symmetric_overlap(
1206            [1.3, 0.0].into(),
1207            &square,
1208            &triangle,
1209            r_square,
1210            r_triangle,
1211            false,
1212        );
1213
1214        assert_symmetric_overlap(
1215            [-1.3, 0.0].into(),
1216            &square,
1217            &triangle,
1218            r_square,
1219            r_triangle,
1220            false,
1221        );
1222
1223        assert_symmetric_overlap(
1224            [0.0, 1.3].into(),
1225            &square,
1226            &triangle,
1227            r_square,
1228            r_triangle,
1229            false,
1230        );
1231
1232        assert_symmetric_overlap(
1233            [0.0, -1.3].into(),
1234            &square,
1235            &triangle,
1236            r_square,
1237            r_triangle,
1238            false,
1239        );
1240
1241        assert_symmetric_overlap(
1242            [1.2, 0.2].into(),
1243            &square,
1244            &triangle,
1245            r_square,
1246            r_triangle,
1247            true,
1248        );
1249
1250        assert_symmetric_overlap(
1251            [-0.7, -0.2].into(),
1252            &square,
1253            &triangle,
1254            r_square,
1255            r_triangle,
1256            true,
1257        );
1258
1259        assert_symmetric_overlap(
1260            [0.4, 1.1].into(),
1261            &square,
1262            &triangle,
1263            r_square,
1264            r_triangle,
1265            true,
1266        );
1267
1268        assert_symmetric_overlap(
1269            [-0.2, -1.2].into(),
1270            &square,
1271            &triangle,
1272            r_square,
1273            r_triangle,
1274            true,
1275        );
1276    }
1277}