1use 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#[derive(Clone, Debug, PartialEq, Serialize, Deserialize)]
81pub struct ConvexSurfaceMesh2d {
82 vertices: Vec<Cartesian<2>>,
84 bounding_radius: PositiveReal,
86}
87
88#[inline]
90fn find_lowest_leftmost(vertices: &[Cartesian<2>]) -> Option<usize> {
91 vertices.iter().position_min_by(|a, b| {
92 a[1].total_cmp(&b[1]).then(a[0].total_cmp(&b[0]))
94 })
95}
96
97#[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#[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 #[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 #[inline]
202 pub fn construct_convex_hull(
203 mut points: Vec<Cartesian<2>>,
204 ) -> Result<Vec<Cartesian<2>>, Error> {
205 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 points.swap(0, anchor_idx);
214 let anchor = points[0];
215
216 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 let mut n_vertices_on_hull = 2;
226 let mut next_candidate = 2;
227
228 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 n_vertices_on_hull -= 1;
238 } else {
239 break;
240 }
241 }
242 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 #[inline]
258 #[must_use]
259 pub fn vertices(&self) -> &[Cartesian<2>] {
260 &self.vertices
261 }
262}
263
264impl SupportMapping<Cartesian<2>> for ConvexSurfaceMesh2d {
265 #[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 #[inline]
352 fn bounding_sphere_radius(&self) -> PositiveReal {
353 self.bounding_radius
354 }
355}
356
357impl Volume for ConvexSurfaceMesh2d {
358 #[inline]
380 fn volume(&self) -> f64 {
381 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 #[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 #[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#[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#[inline]
544fn is_separating(a: &ConvexSurfaceMesh2d, p: &Cartesian<2>, n: &Cartesian<2>) -> bool {
545 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 #[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], [1.0, 0.5],
679 [1.0, 1.0], [0.5, 1.0],
681 [0.0, 1.0], [0.0, 0.5], ]
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 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 for i in 0..20 {
705 pts.push([f64::from(i) / 19.0, 0.0]);
706 }
707 for i in 0..20 {
709 pts.push([1.0, f64::from(i) / 19.0]);
710 }
711 for i in 0..20 {
713 pts.push([f64::from(i) / 19.0, 1.0]);
714 }
715 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 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 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 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 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 assert!(vertices.len() >= 3);
791 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], [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 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], [0.75, 1.0], ]
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 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], [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 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 mut points: Vec<[f64; 2]> = (0..10).map(|i| [f64::from(i) * 2.0 / 9.0, 0.0]).collect();
926 points.push([1.0, 1.0]); let points: Vec<Cartesian<2>> = points.into_iter().map(Cartesian::from).collect();
928 let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
929 .expect("hard-coded points should lie on a convex hull");
930 assert!(vertices.len() >= 3);
932 }
933
934 #[rstest]
935 fn test_collinear_from_anchor() {
936 let points: Vec<Cartesian<2>> = vec![
937 [0.0, 0.0], [1.0, 1.0],
939 [2.0, 2.0],
940 [3.0, 3.0], [1.0, 0.0],
942 [0.0, 1.0], ]
944 .into_iter()
945 .map(Cartesian::from)
946 .collect();
947 let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
948 .expect("hard-coded points should lie on a convex hull");
949 assert!(vertices.len() >= 3);
951
952 let hull = [
953 [0.0, 0.0].into(),
954 [1.0, 0.0].into(),
955 [3.0, 3.0].into(),
956 [0.0, 1.0].into(),
957 ];
958 itertools::assert_equal(&vertices, &hull);
959 }
960
961 #[rstest]
962 fn test_radial_collinear_multiple_directions() {
963 let points: Vec<Cartesian<2>> = vec![
964 [0.0, 0.0], [1.0, 0.0],
966 [2.0, 0.0],
967 [3.0, 0.0], [0.0, 1.0],
969 [0.0, 2.0],
970 [0.0, 3.0], [1.0, 1.0],
972 [2.0, 2.0], [-1.0, 0.0],
974 [-2.0, 0.0], ]
976 .into_iter()
977 .map(Cartesian::from)
978 .collect();
979 let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
980 .expect("hard-coded points should lie on a convex hull");
981 assert!(vertices.len() >= 3);
983 }
984
985 #[rstest]
986 fn test_star_pattern() {
987 let mut points: Vec<Cartesian<2>> = vec![Cartesian::from([0.0, 0.0])]; for i in 0..8 {
990 let angle = 2.0 * std::f64::consts::PI * f64::from(i) / 8.0;
991 for r in [0.5, 1.0, 1.5] {
992 points.push(Cartesian::from([r * angle.cos(), r * angle.sin()]));
993 }
994 }
995 let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
996 .expect("hard-coded points should lie on a convex hull");
997 assert!(vertices.len() >= 8); }
1000
1001 #[rstest]
1002 fn test_minimum_triangle() {
1003 let points: Vec<Cartesian<2>> = vec![[0.0, 0.0], [1.0, 0.0], [0.5, 1.0]]
1004 .into_iter()
1005 .map(Cartesian::from)
1006 .collect();
1007 let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
1008 .expect("hard-coded points should lie on a convex hull");
1009 assert_eq!(vertices.len(), 3);
1010
1011 let hull = [[0.0, 0.0].into(), [1.0, 0.0].into(), [0.5, 1.0].into()];
1012 itertools::assert_equal(&vertices, &hull);
1013 }
1014
1015 #[rstest]
1016 fn test_negative_coordinates() {
1017 let points: Vec<Cartesian<2>> = vec![[1.0, 1.0], [1.0, -1.0], [-1.0, 1.0], [-1.0, -1.0]]
1018 .into_iter()
1019 .map(Cartesian::from)
1020 .collect();
1021 let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
1022 .expect("hard-coded points should lie on a convex hull");
1023 assert_eq!(vertices.len(), 4);
1024
1025 let hull = [
1026 [-1.0, -1.0].into(),
1027 [1.0, -1.0].into(),
1028 [1.0, 1.0].into(),
1029 [-1.0, 1.0].into(),
1030 ];
1031 itertools::assert_equal(&vertices, &hull);
1032 }
1033
1034 #[rstest]
1035 fn test_large_coordinates() {
1036 let points: Vec<Cartesian<2>> = vec![[0.0, 0.0], [1e6, 0.0], [1e6, 1e6], [0.0, 1e6]]
1037 .into_iter()
1038 .map(Cartesian::from)
1039 .collect();
1040 let vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)
1041 .expect("hard-coded points should lie on a convex hull");
1042 assert_eq!(vertices.len(), 4);
1043
1044 let hull = [
1045 [0.0, 0.0].into(),
1046 [1e6, 0.0].into(),
1047 [1e6, 1e6].into(),
1048 [0.0, 1e6].into(),
1049 ];
1050 itertools::assert_equal(&vertices, &hull);
1051 }
1052
1053 #[rstest]
1054 fn test_degenerate() {
1055 let points: Vec<Cartesian<2>> = vec![[0.0, 0.0], [0.5, 0.5], [0.25, 0.25], [1.0, 1.0]]
1056 .into_iter()
1057 .map(Cartesian::from)
1058 .collect();
1059 let result = ConvexSurfaceMesh2d::construct_convex_hull(points);
1060 check!(result == Err(Error::DegeneratePolytope));
1061 }
1062
1063 #[test]
1064 fn support_mapping_2d() {
1065 let cuboid = ConvexSurfaceMesh2d::from_point_set([
1066 [-1.0, -2.0].into(),
1067 [1.0, -2.0].into(),
1068 [1.0, 2.0].into(),
1069 [-1.0, 2.0].into(),
1070 ])
1071 .expect("hard-coded vertices form a polygon");
1072
1073 assert_relative_eq!(
1074 cuboid.support_mapping(&Cartesian::from([1.0, 0.1])),
1075 [1.0, 2.0].into()
1076 );
1077 assert_relative_eq!(
1078 cuboid.support_mapping(&Cartesian::from([1.0, -0.1])),
1079 [1.0, -2.0].into()
1080 );
1081 assert_relative_eq!(
1082 cuboid.support_mapping(&Cartesian::from([-0.1, 1.0])),
1083 [-1.0, 2.0].into()
1084 );
1085 assert_relative_eq!(
1086 cuboid.support_mapping(&Cartesian::from([-0.1, -1.0])),
1087 [-1.0, -2.0].into()
1088 );
1089 }
1090
1091 #[fixture]
1094 fn square() -> ConvexSurfaceMesh2d {
1095 ConvexSurfaceMesh2d::from_point_set([
1096 [-0.5, -0.5].into(),
1097 [0.5, -0.5].into(),
1098 [0.5, 0.5].into(),
1099 [-0.5, 0.5].into(),
1100 ])
1101 .expect("hard-coded vertices form a valid polygon")
1102 }
1103
1104 #[fixture]
1105 fn triangle() -> ConvexSurfaceMesh2d {
1106 ConvexSurfaceMesh2d::from_point_set([
1107 [-0.5, -0.5].into(),
1108 [0.5, -0.5].into(),
1109 [0.5, 0.5].into(),
1110 ])
1111 .expect("hard-coded vertices form a valid polygon")
1112 }
1113
1114 #[rstest]
1115 fn square_no_rot(square: ConvexSurfaceMesh2d) {
1116 let a = Angle::identity();
1117 assert!(!square.intersects_at(&square, &[10.0, 0.0].into(), &a));
1118 assert!(!square.intersects_at(&square, &[-10.0, 0.0].into(), &a));
1119
1120 assert!(!square.intersects_at(&square, &[1.1, 0.0].into(), &a));
1121 assert!(!square.intersects_at(&square, &[-1.1, 0.0].into(), &a));
1122 assert!(!square.intersects_at(&square, &[0.0, 1.1].into(), &a));
1123 assert!(!square.intersects_at(&square, &[0.0, -1.1].into(), &a));
1124
1125 assert!(square.intersects_at(&square, &[0.9, 0.2].into(), &a));
1126 assert!(square.intersects_at(&square, &[-0.9, 0.2].into(), &a));
1127 assert!(square.intersects_at(&square, &[-0.2, 0.9].into(), &a));
1128 assert!(square.intersects_at(&square, &[-0.2, -0.9].into(), &a));
1129
1130 assert!(square.intersects_at(&square, &[1.0, 0.2].into(), &a));
1131 }
1132
1133 #[rstest]
1134 fn square_rot(square: ConvexSurfaceMesh2d) {
1135 let a = Angle::from(PI / 4.0);
1136
1137 assert!(!square.intersects_at(&square, &[10.0, 0.0].into(), &a));
1138 assert!(!square.intersects_at(&square, &[-10.0, 0.0].into(), &a));
1139
1140 assert!(!square.intersects_at(&square, &[1.3, 0.0].into(), &a));
1141 assert!(!square.intersects_at(&square, &[-1.3, 0.0].into(), &a));
1142 assert!(!square.intersects_at(&square, &[0.0, 1.3].into(), &a));
1143 assert!(!square.intersects_at(&square, &[0.0, -1.3].into(), &a));
1144
1145 assert!(!square.intersects_at(&square, &[1.3, 0.2].into(), &a));
1146 assert!(!square.intersects_at(&square, &[-1.3, 0.2].into(), &a));
1147 assert!(!square.intersects_at(&square, &[-0.2, 1.3].into(), &a));
1148 assert!(!square.intersects_at(&square, &[-0.2, -1.3].into(), &a));
1149
1150 assert!(square.intersects_at(&square, &[1.2, 0.2].into(), &a));
1151 assert!(square.intersects_at(&square, &[-1.2, 0.2].into(), &a));
1152 assert!(square.intersects_at(&square, &[-0.2, 1.2].into(), &a));
1153 assert!(square.intersects_at(&square, &[-0.2, -1.2].into(), &a));
1154 }
1155
1156 fn test_overlap<A, B, R, const N: usize>(
1157 r_ab: Cartesian<N>,
1158 a: &A,
1159 b: &B,
1160 o_a: R,
1161 o_b: &R,
1162 ) -> bool
1163 where
1164 R: Rotation + Rotate<Cartesian<N>>,
1165 A: IntersectsAt<B, Cartesian<N>, R>,
1166 {
1167 let r_a_inverted = o_a.inverted();
1168 let v_ij = r_a_inverted.rotate(&r_ab);
1169 let o_ij = o_b.combine(&r_a_inverted);
1170 a.intersects_at(b, &v_ij, &o_ij)
1171 }
1172
1173 fn assert_symmetric_overlap<A, B, R, const N: usize>(
1174 r_ab: Cartesian<N>,
1175 a: &A,
1176 b: &B,
1177 r_a: R,
1178 r_b: R,
1179 expected: bool,
1180 ) where
1181 R: Rotation + Rotate<Cartesian<N>>,
1182 A: IntersectsAt<B, Cartesian<N>, R>,
1183 B: IntersectsAt<A, Cartesian<N>, R>,
1184 {
1185 assert_eq!(test_overlap(r_ab, a, b, r_a, &r_b), expected);
1186 assert_eq!(test_overlap(-r_ab, b, a, r_b, &r_a), expected);
1187 }
1188
1189 #[rstest]
1190 fn square_triangle(square: ConvexSurfaceMesh2d, triangle: ConvexSurfaceMesh2d) {
1191 let r_square = Angle::from(-PI / 4.0);
1192 let r_triangle = Angle::from(PI);
1193
1194 assert_symmetric_overlap(
1195 [10.0, 0.0].into(),
1196 &square,
1197 &triangle,
1198 r_square,
1199 r_triangle,
1200 false,
1201 );
1202
1203 assert_symmetric_overlap(
1204 [1.3, 0.0].into(),
1205 &square,
1206 &triangle,
1207 r_square,
1208 r_triangle,
1209 false,
1210 );
1211
1212 assert_symmetric_overlap(
1213 [-1.3, 0.0].into(),
1214 &square,
1215 &triangle,
1216 r_square,
1217 r_triangle,
1218 false,
1219 );
1220
1221 assert_symmetric_overlap(
1222 [0.0, 1.3].into(),
1223 &square,
1224 &triangle,
1225 r_square,
1226 r_triangle,
1227 false,
1228 );
1229
1230 assert_symmetric_overlap(
1231 [0.0, -1.3].into(),
1232 &square,
1233 &triangle,
1234 r_square,
1235 r_triangle,
1236 false,
1237 );
1238
1239 assert_symmetric_overlap(
1240 [1.2, 0.2].into(),
1241 &square,
1242 &triangle,
1243 r_square,
1244 r_triangle,
1245 true,
1246 );
1247
1248 assert_symmetric_overlap(
1249 [-0.7, -0.2].into(),
1250 &square,
1251 &triangle,
1252 r_square,
1253 r_triangle,
1254 true,
1255 );
1256
1257 assert_symmetric_overlap(
1258 [0.4, 1.1].into(),
1259 &square,
1260 &triangle,
1261 r_square,
1262 r_triangle,
1263 true,
1264 );
1265
1266 assert_symmetric_overlap(
1267 [-0.2, -1.2].into(),
1268 &square,
1269 &triangle,
1270 r_square,
1271 r_triangle,
1272 true,
1273 );
1274 }
1275}