use robust_geo as rg;
pub use nalgebra;
use nalgebra::{Vector1, Vector2, Vector3};
type Vec1 = Vector1<f64>;
type Vec2 = Vector2<f64>;
type Vec3 = Vector3<f64>;
macro_rules! sorted_fn {
($name:ident, $n:expr) => {
fn $name<Idx: Ord + Copy>(mut arr: [Idx; $n]) -> ([Idx; $n], bool) {
let mut num_swaps = 0;
for i in 1..$n {
for j in (0..i).rev() {
if arr[j] > arr[j + 1] {
arr.swap(j, j + 1);
num_swaps += 1;
} else {
break;
}
}
}
(arr, num_swaps % 2 != 0)
}
};
}
sorted_fn!(sorted_3, 3);
sorted_fn!(sorted_4, 4);
sorted_fn!(sorted_5, 5);
pub fn orient_1d<T: ?Sized, Idx: Ord + Copy>(
list: &T,
index_fn: impl Fn(&T, Idx) -> Vec1,
i: Idx,
j: Idx,
) -> bool {
let pi = index_fn(list, i);
let pj = index_fn(list, j);
pi > pj || (pi == pj && i < j)
}
macro_rules! case {
(2: $pi:ident, $pj:ident, @ m2, != $odd:expr) => {
let val = rg::magnitude_cmp_2d($pi, $pj);
if val != 0.0 {
return (val > 0.0) != $odd;
}
};
(2: $pi:ident, $pj:ident, @ m3, != $odd:expr) => {
let val = rg::magnitude_cmp_3d($pi, $pj);
if val != 0.0 {
return (val > 0.0) != $odd;
}
};
(2: $pi:ident, $pj:ident, $(@ $swiz:ident,)? != $odd:expr) => {
if $pi$(.$swiz)? != $pj$(.$swiz)? {
return ($pi$(.$swiz)? > $pj$(.$swiz)?) != $odd;
}
};
(3: $pi:ident, $pj:ident, $pk:ident, @ $swiz:ident m2, != $odd:expr) => {
let val = rg::sign_det_x_x2y2($pi.$swiz(), $pj.$swiz(), $pk.$swiz());
if val != 0.0 {
return (val > 0.0) != $odd;
}
};
(3: $pi:ident, $pj:ident, $pk:ident, @ $swiz:ident m3, != $odd:expr) => {
let val = rg::sign_det_x_x2y2z2($pi.$swiz(), $pj.$swiz(), $pk.$swiz());
if val != 0.0 {
return (val > 0.0) != $odd;
}
};
(3: $pi:ident, $pj:ident, $pk:ident, $(@ $swiz:ident,)? != $odd:expr) => {
let val = rg::orient_2d($pi$(.$swiz())?, $pj$(.$swiz())?, $pk$(.$swiz())?);
if val != 0.0 {
return (val > 0.0) != $odd;
}
};
(4: $pi:ident, $pj:ident, $pk:ident, $pl:ident, @ xy m2, != $odd:expr) => {
let val = rg::in_circle($pi, $pj, $pk, $pl);
if val != 0.0 {
return (val > 0.0) != $odd;
}
};
(4: $pi:ident, $pj:ident, $pk:ident, $pl:ident, @ $swiz:ident m3, != $odd:expr) => {
let val = rg::sign_det_x_y_x2y2z2($pi.$swiz(), $pj.$swiz(), $pk.$swiz(), $pl.$swiz());
if val != 0.0 {
return (val > 0.0) != $odd;
}
};
(4: $pi:ident, $pj:ident, $pk:ident, $pl:ident, $(@ $swiz:ident,)? != $odd:expr) => {
let val = rg::orient_3d($pi$(.$swiz())?, $pj$(.$swiz())?, $pk$(.$swiz())?, $pl$(.$swiz())?);
if val != 0.0 {
return (val > 0.0) != $odd;
}
};
(5: $pi:ident, $pj:ident, $pk:ident, $pl:ident, $pm:ident, @ xyz m3, != $odd:expr) => {
let val = rg::in_sphere($pi, $pj, $pk, $pl, $pm);
if val != 0.0 {
return (val > 0.0) != $odd;
}
};
}
pub fn orient_2d<T: ?Sized, Idx: Ord + Copy>(
list: &T,
index_fn: impl Fn(&T, Idx) -> Vec2,
i: Idx,
j: Idx,
k: Idx,
) -> bool {
let ([i, j, k], odd) = sorted_3([i, j, k]);
let pi = index_fn(list, i);
let pj = index_fn(list, j);
let pk = index_fn(list, k);
case!(3: pi, pj, pk, != odd);
case!(2: pk, pj, @ x, != odd);
case!(2: pj, pk, @ y, != odd);
case!(2: pi, pk, @ x, != odd);
!odd
}
pub fn orient_3d<T: ?Sized, Idx: Ord + Copy>(
list: &T,
index_fn: impl Fn(&T, Idx) -> Vec3,
i: Idx,
j: Idx,
k: Idx,
l: Idx,
) -> bool {
let ([i, j, k, l], odd) = sorted_4([i, j, k, l]);
let pi = index_fn(list, i);
let pj = index_fn(list, j);
let pk = index_fn(list, k);
let pl = index_fn(list, l);
case!(4: pi, pj, pk, pl, != odd);
case!(3: pj, pk, pl, @ xy, != odd);
case!(3: pj, pk, pl, @ zx, != odd);
case!(3: pj, pk, pl, @ yz, != odd);
case!(3: pi, pk, pl, @ yx, != odd);
case!(2: pk, pl, @ x, != odd);
case!(2: pl, pk, @ y, != odd);
case!(3: pi, pk, pl, @ xz, != odd);
case!(2: pk, pl, @ z, != odd);
case!(3: pi, pj, pl, @ xy, != odd);
case!(2: pl, pj, @ x, != odd);
case!(2: pj, pl, @ y, != odd);
case!(2: pi, pl, @ x, != odd);
!odd
}
pub fn in_circle<T: ?Sized, Idx: Ord + Copy>(
list: &T,
index_fn: impl Fn(&T, Idx) -> Vec2 + Clone,
i: Idx,
j: Idx,
k: Idx,
l: Idx,
) -> bool {
simplicity_derive::generate_in_hypersphere!{list, index_fn, i, j, k, l}
}
pub fn in_circle_unoriented<T: ?Sized, Idx: Ord + Copy>(
list: &T,
index_fn: impl Fn(&T, Idx) -> Vec2 + Clone,
i: Idx,
j: Idx,
k: Idx,
l: Idx,
) -> bool {
orient_2d(list, index_fn.clone(), i, j, k) == in_circle(list, index_fn, i, j, k, l)
}
pub fn in_sphere<T: ?Sized, Idx: Ord + Copy>(
list: &T,
index_fn: impl Fn(&T, Idx) -> Vec3 + Clone,
i: Idx,
j: Idx,
k: Idx,
l: Idx,
m: Idx,
) -> bool {
simplicity_derive::generate_in_hypersphere!{list, index_fn, i, j, k, l, m}
}
pub fn in_sphere_unoriented<T: ?Sized, Idx: Ord + Copy>(
list: &T,
index_fn: impl Fn(&T, Idx) -> Vec3 + Clone,
i: Idx,
j: Idx,
k: Idx,
l: Idx,
m: Idx,
) -> bool {
orient_3d(list, index_fn.clone(), i, j, k, l) == in_sphere(list, index_fn, i, j, k, l, m)
}
#[cfg(test)]
mod tests {
use super::*;
use test_case::test_case;
macro_rules! case {
($arr:expr => $pi:ident, $pj:ident, @ m2) => {
let val = rg::magnitude_cmp_2d($pi, $pj);
if val != 0.0 {
return $arr;
}
};
($arr:expr => $pi:ident, $pj:ident, @ m3) => {
let val = rg::magnitude_cmp_3d($pi, $pj);
if val != 0.0 {
return $arr;
}
};
($arr:expr => $pi:ident, $pj:ident $(, @ $swiz:ident)?) => {
if $pi$(.$swiz)? != $pj$(.$swiz)? {
return $arr;
}
};
($arr:expr => $pi:ident, $pj:ident, $pk:ident, @ $swiz:ident m2) => {
let val = rg::sign_det_x_x2y2($pi.$swiz(), $pj.$swiz(), $pk.$swiz());
if val != 0.0 {
return $arr;
}
};
($arr:expr => $pi:ident, $pj:ident, $pk:ident, @ $swiz:ident m3) => {
let val = rg::sign_det_x_x2y2z2($pi.$swiz(), $pj.$swiz(), $pk.$swiz());
if val != 0.0 {
return $arr;
}
};
($arr:expr => $pi:ident, $pj:ident, $pk:ident $(, @ $swiz:ident)?) => {
let val = rg::orient_2d($pi$(.$swiz())?, $pj$(.$swiz())?, $pk$(.$swiz())?);
if val != 0.0 {
return $arr;
}
};
($arr:expr => $pi:ident, $pj:ident, $pk:ident, $pl:ident, @ xy m2) => {
let val = rg::in_circle($pi, $pj, $pk, $pl);
if val != 0.0 {
return $arr;
}
};
($arr:expr => $pi:ident, $pj:ident, $pk:ident, $pl:ident, @ $swiz:ident m3) => {
let val = rg::sign_det_x_y_x2y2z2($pi.$swiz(), $pj.$swiz(), $pk.$swiz(), $pl.$swiz());
if val != 0.0 {
return $arr;
}
};
($arr:expr => $pi:ident, $pj:ident, $pk:ident, $pl:ident $(, @ $swiz:ident)?) => {
let val = rg::orient_3d($pi$(.$swiz())?, $pj$(.$swiz())?, $pk$(.$swiz())?, $pl$(.$swiz())?);
if val != 0.0 {
return $arr;
}
};
($arr:expr => $pi:ident, $pj:ident, $pk:ident, $pl:ident, $pm:ident, @ xyz m3) => {
let val = rg::in_sphere($pi, $pj, $pk, $pl, $pm);
if val != 0.0 {
return $arr;
}
};
}
pub fn orient_2d_case<T: ?Sized>(
list: &T,
index_fn: impl Fn(&T, usize) -> Vec2,
i: usize,
j: usize,
k: usize,
) -> [usize; 3] {
let ([i, j, k], _) = sorted_3([i, j, k]);
let pi = index_fn(list, i);
let pj = index_fn(list, j);
let pk = index_fn(list, k);
case!([3, 3, 3] => pi, pj, pk);
case!([2, 3, 3] => pk, pj, @ x);
case!([1, 3, 3] => pj, pk, @ y);
case!([2, 2, 3] => pi, pk, @ x);
[1, 2, 3]
}
pub fn orient_3d_case<T: ?Sized>(
list: &T,
index_fn: impl Fn(&T, usize) -> Vec3,
i: usize,
j: usize,
k: usize,
l: usize,
) -> [usize; 4] {
let ([i, j, k, l], _) = sorted_4([i, j, k, l]);
let pi = index_fn(list, i);
let pj = index_fn(list, j);
let pk = index_fn(list, k);
let pl = index_fn(list, l);
case!([4, 4, 4, 4] => pi, pj, pk, pl);
case!([3, 4, 4, 4] => pj, pk, pl, @ xy);
case!([2, 4, 4, 4] => pj, pk, pl, @ zx);
case!([1, 4, 4, 4] => pj, pk, pl, @ yz);
case!([3, 3, 4, 4] => pi, pk, pl, @ yx);
case!([2, 3, 4, 4] => pk, pl, @ x);
case!([1, 3, 4, 4] => pl, pk, @ y);
case!([2, 2, 4, 4] => pi, pk, pl, @ xz);
case!([1, 2, 4, 4] => pk, pl, @ z);
case!([3, 3, 3, 4] => pi, pj, pl, @ xy);
case!([2, 3, 3, 4] => pl, pj, @ x);
case!([1, 3, 3, 4] => pj, pl, @ y);
case!([2, 2, 3, 4] => pi, pl, @ x);
[1, 2, 3, 4]
}
#[test]
fn orient_1d_positive() {
let points = vec![0.0, 1.0];
assert!(orient_1d(&points, |l, i| Vector1::new(l[i]), 1, 0))
}
#[test]
fn orient_1d_negative() {
let points = vec![0.0, 1.0];
assert!(!orient_1d(&points, |l, i| Vector1::new(l[i]), 0, 1))
}
#[test]
fn orient_1d_positive_degenerate() {
let points = vec![0.0, 0.0];
assert!(orient_1d(&points, |l, i| Vector1::new(l[i]), 0, 1))
}
#[test]
fn orient_1d_negative_degenerate() {
let points = vec![0.0, 0.0];
assert!(!orient_1d(&points, |l, i| Vector1::new(l[i]), 1, 0))
}
#[test_case([[0.0, 0.0], [1.0, 0.0], [2.0, 1.0]], [3,3,3] ; "General")]
#[test_case([[0.0, 0.0], [1.0, 1.0], [2.0, 2.0]], [2,3,3] ; "Collinear")]
#[test_case([[0.0, 0.0], [0.0, 2.0], [0.0, 1.0]], [1,3,3] ; "Collinear, pj.x = pk.x")]
#[test_case([[1.0, 0.0], [0.0, 2.0], [0.0, 2.0]], [2,2,3] ; "pj = pk")]
#[test_case([[0.0, 0.0], [0.0, 2.0], [0.0, 2.0]], [1,2,3] ; "pj = pk, pi.x = pk.x")]
fn test_orient_2d(points: [[f64; 2]; 3], case: [usize; 3]) {
let points = points
.iter()
.copied()
.map(Vector2::from)
.collect::<Vec<_>>();
assert!(orient_2d(&points, |l, i| l[i], 0, 1, 2));
assert!(!orient_2d(&points, |l, i| l[i], 0, 2, 1));
assert!(!orient_2d(&points, |l, i| l[i], 1, 0, 2));
assert!(orient_2d(&points, |l, i| l[i], 1, 2, 0));
assert!(orient_2d(&points, |l, i| l[i], 2, 0, 1));
assert!(!orient_2d(&points, |l, i| l[i], 2, 1, 0));
assert_eq!(orient_2d_case(&points, |l, i| l[i], 0, 1, 2), case);
}
#[test_case([[0.0, 0.0, 0.0], [0.0, 0.0, 1.0], [0.0, 1.0, 0.0], [1.0, 0.0, 0.0]], [4,4,4,4] ; "General")]
#[test_case([[0.0, 0.0, 0.0], [1.0, 1.0, 1.0], [3.0, 4.0, 5.0], [2.0, 3.0, 4.0]], [3,4,4,4] ; "Coplanar")]
#[test_case([[0.0, 0.0, 0.0], [1.0, 1.0, 1.0], [2.0, 2.0, 4.0], [3.0, 3.0, 5.0]], [2,4,4,4] ; "Coplanar, pj pk pl @ xy collinear")]
#[test_case([[1.0, 0.0, 0.0], [1.0, 1.0, 1.0], [1.0, 4.0, 2.0], [1.0, 5.0, 3.0]], [1,4,4,4] ; "Coplanar, pj.x = pk.x = pl.x or pj pk pl collinear")]
#[test_case([[0.0, 0.0, 0.0], [1.0, 2.0, 3.0], [2.0, 3.0, 4.0], [3.0, 4.0, 5.0]], [3,3,4,4] ; "pj pk pl collinear")]
#[test_case([[0.0, 0.0, 0.0], [1.0, 1.0, 3.0], [3.0, 3.0, 5.0], [2.0, 2.0, 4.0]], [2,3,4,4] ; "pj pk pl collinear, pi pk pl @ xy collinear")]
#[test_case([[0.0, 0.0, 0.0], [0.0, 1.0, 3.0], [0.0, 2.0, 4.0], [0.0, 3.0, 5.0]], [1,3,4,4] ; "pj pk pl collinear, pi pk pl @ xy collinear, pk.x = pl.x")]
#[test_case([[1.0, 0.0, 0.0], [0.0, 2.0, 3.0], [0.0, 2.0, 5.0], [0.0, 2.0, 4.0]], [2,2,4,4] ; "pj pk pl collinear, pi pk pl @ xy collinear, pk.xy = pl.xy")]
#[test_case([[0.0, 0.0, 0.0], [0.0, 2.0, 3.0], [0.0, 2.0, 4.0], [0.0, 2.0, 3.0]], [1,2,4,4] ; "pj pk pl collinear, pi.x = pk.x = pl.x or pi pk pl collinear, pk.xy = pl.xy")]
#[test_case([[0.0, 0.0, 0.0], [1.0, 0.0, 0.0], [2.0, 1.0, 0.0], [2.0, 1.0, 0.0]], [3,3,3,4] ; "pk = pl")]
#[test_case([[0.0, 0.0, 0.0], [1.0, 1.0, 0.0], [2.0, 2.0, 0.0], [2.0, 2.0, 0.0]], [2,3,3,4] ; "pk = pl, pi pj pk @ xy collinear")]
#[test_case([[0.0, 0.0, 0.0], [0.0, 2.0, 0.0], [0.0, 1.0, 0.0], [0.0, 1.0, 0.0]], [1,3,3,4] ; "pk = pl, pi pj pk @ xy collinear, pj.x = pk.x")]
#[test_case([[1.0, 0.0, 0.0], [0.0, 2.0, 0.0], [0.0, 2.0, 0.0], [0.0, 2.0, 0.0]], [2,2,3,4] ; "pk = pl, pi pj pk @ xy collinear, pj.xy = pk.xy")]
#[test_case([[0.0, 0.0, 0.0], [0.0, 2.0, 0.0], [0.0, 2.0, 0.0], [0.0, 2.0, 0.0]], [1,2,3,4] ; "pk = pl, pi pj pk @ xy collinear, pj.xy = pk.xy, pi.x = pk.x")]
fn test_orient_3d(points: [[f64; 3]; 4], case: [usize; 4]) {
let points = points
.iter()
.copied()
.map(Vector3::from)
.collect::<Vec<_>>();
assert!(orient_3d(&points, |l, i| l[i], 0, 1, 2, 3));
assert!(!orient_3d(&points, |l, i| l[i], 3, 2, 0, 1));
assert_eq!(orient_3d_case(&points, |l, i| l[i], 0, 1, 2, 3), case);
}
#[test]
fn test_in_circle_unoriented_general() {
let points = [[0.0, 0.0], [0.0, 2.0], [2.0, 2.0], [1.0, 1.0]];
let points = points
.iter()
.copied()
.map(Vector2::from)
.collect::<Vec<_>>();
assert!(in_circle_unoriented(&points, |l, i| l[i], 0, 1, 2, 3));
assert!(in_circle_unoriented(&points, |l, i| l[i], 0, 2, 1, 3));
assert!(in_circle_unoriented(&points, |l, i| l[i], 1, 2, 0, 3));
assert!(in_circle_unoriented(&points, |l, i| l[i], 1, 0, 2, 3));
assert!(in_circle_unoriented(&points, |l, i| l[i], 2, 0, 1, 3));
assert!(in_circle_unoriented(&points, |l, i| l[i], 2, 1, 0, 3));
assert!(
(in_circle_unoriented(&points, |l, i| l[i], 0, 1, 2, 3)
== in_circle_unoriented(&points, |l, i| l[i], 0, 1, 3, 2))
== (orient_2d(&points, |l, i| l[i], 0, 1, 3)
!= orient_2d(&points, |l, i| l[i], 0, 1, 2))
);
}
#[test]
fn test_in_circle_unoriented_cocircular() {
let points = [[0.0, 0.0], [0.0, 0.0], [1.0, 0.0], [0.0, 1.0]];
let points = points
.iter()
.copied()
.map(Vector2::from)
.collect::<Vec<_>>();
assert!(in_circle_unoriented(&points, |l, i| l[i], 1, 2, 3, 0));
assert!(in_circle_unoriented(&points, |l, i| l[i], 1, 3, 2, 0));
assert!(in_circle_unoriented(&points, |l, i| l[i], 2, 3, 1, 0));
assert!(in_circle_unoriented(&points, |l, i| l[i], 1, 2, 3, 0));
assert!(in_circle_unoriented(&points, |l, i| l[i], 3, 1, 2, 0));
assert!(in_circle_unoriented(&points, |l, i| l[i], 3, 2, 1, 0));
assert!(
(in_circle_unoriented(&points, |l, i| l[i], 0, 1, 2, 3)
== in_circle_unoriented(&points, |l, i| l[i], 0, 1, 3, 2))
== (orient_2d(&points, |l, i| l[i], 0, 1, 3)
!= orient_2d(&points, |l, i| l[i], 0, 1, 2))
);
}
#[test]
fn test_in_sphere_unoriented_general() {
let points = [[0,0,0], [4,0,0], [0,4,0], [0,0,4], [1,1,1]];
let points = points
.iter()
.copied()
.map(|[x, y, z]| Vector3::new(x as f64, y as f64, z as f64))
.collect::<Vec<_>>();
assert!(in_sphere_unoriented(&points, |l, i| l[i], 0, 1, 2, 3, 4));
assert!(in_sphere_unoriented(&points, |l, i| l[i], 0, 2, 1, 3, 4));
assert!(in_sphere_unoriented(&points, |l, i| l[i], 1, 2, 0, 3, 4));
assert!(in_sphere_unoriented(&points, |l, i| l[i], 1, 3, 0, 2, 4));
assert!(in_sphere_unoriented(&points, |l, i| l[i], 2, 3, 0, 1, 4));
assert!(in_sphere_unoriented(&points, |l, i| l[i], 2, 3, 1, 0, 4));
assert!(
(in_sphere_unoriented(&points, |l, i| l[i], 0, 1, 2, 3, 4)
== in_sphere_unoriented(&points, |l, i| l[i], 0, 1, 2, 4, 3))
== (orient_3d(&points, |l, i| l[i], 0, 1, 2, 3)
!= orient_3d(&points, |l, i| l[i], 0, 1, 2, 4))
);
}
#[test]
fn test_in_sphere_unoriented_cospherical() {
let points = [[0,0,0], [0,0,0], [1,0,0], [0,0,1], [0,1,0]];
let points = points
.iter()
.copied()
.map(|[x, y, z]| Vector3::new(x as f64, y as f64, z as f64))
.collect::<Vec<_>>();
assert!(in_sphere_unoriented(&points, |l, i| l[i], 1, 2, 3, 4, 0));
assert!(in_sphere_unoriented(&points, |l, i| l[i], 1, 3, 2, 4, 0));
assert!(in_sphere_unoriented(&points, |l, i| l[i], 2, 3, 1, 4, 0));
assert!(in_sphere_unoriented(&points, |l, i| l[i], 2, 4, 1, 3, 0));
assert!(in_sphere_unoriented(&points, |l, i| l[i], 3, 4, 1, 2, 0));
assert!(in_sphere_unoriented(&points, |l, i| l[i], 3, 4, 2, 1, 0));
assert!(
(in_sphere_unoriented(&points, |l, i| l[i], 0, 1, 2, 3, 4)
== in_sphere_unoriented(&points, |l, i| l[i], 0, 1, 2, 4, 3))
== (orient_3d(&points, |l, i| l[i], 0, 1, 2, 3)
!= orient_3d(&points, |l, i| l[i], 0, 1, 2, 4))
);
}
}