use std::f64::consts::PI;
use nalgebra::Matrix2;
use crate::geometry::diagram::RegionMask;
use crate::geometry::primitives::Point;
use crate::geometry::projective::Conic;
use crate::geometry::shapes::{Circle, Polygon, Rectangle};
use crate::geometry::traits::{
Area, BoundingBox, Centroid, Closed, DiagramShape, Perimeter, Polygonize,
};
#[derive(Debug, Clone, Copy, PartialEq)]
pub struct Ellipse {
center: Point,
semi_major: f64,
semi_minor: f64,
rotation: f64, }
impl Ellipse {
pub fn new(center: Point, semi_major: f64, semi_minor: f64, rotation: f64) -> Self {
debug_assert!(semi_major > 0.0, "Semi-major axis must be > 0");
debug_assert!(semi_minor > 0.0, "Semi-minor axis must be > 0");
Self {
center,
semi_major,
semi_minor,
rotation,
}
}
pub fn from_radius_ratio(center: Point, radius: f64, log_aspect: f64, rotation: f64) -> Self {
let radius = radius.abs();
let log_aspect = log_aspect.clamp(-1.6094379124341003, 0.0); let aspect_ratio = log_aspect.exp();
let rotation = rotation.rem_euclid(PI);
let semi_major = radius / aspect_ratio.sqrt();
let semi_minor = radius * aspect_ratio.sqrt();
Self::new(center, semi_major, semi_minor, rotation)
}
pub fn from_conic(conic: Conic) -> Option<Self> {
let c = conic.matrix();
let m1 = c[(0, 0)];
let m2 = 2.0 * c[(0, 1)];
let m3 = c[(1, 1)];
let m4 = 2.0 * c[(0, 2)];
let m5 = 2.0 * c[(1, 2)];
let m6 = c[(2, 2)];
if m2 * m2 - 4.0 * m1 * m3 >= 0.0 {
return None;
}
let m = Matrix2::new(m1, m2 / 2.0, m2 / 2.0, m3);
let rhs = nalgebra::Vector2::new(-m4 / 2.0, -m5 / 2.0);
let center_vec = m.lu().solve(&rhs)?;
let xc = center_vec[0];
let yc = center_vec[1];
let f_bar = m6 + m1 * xc * xc + m2 * xc * yc + m3 * yc * yc + m4 * xc + m5 * yc;
if f_bar >= 0.0 {
return None; }
let eig = m.symmetric_eigen();
let lambda1 = eig.eigenvalues[0];
let lambda2 = eig.eigenvalues[1];
let v = eig.eigenvectors;
let (lambda_major, lambda_minor, vec_major) = if lambda1 < lambda2 {
(lambda1, lambda2, v.column(0))
} else {
(lambda2, lambda1, v.column(1))
};
let a2 = -f_bar / lambda_major;
let b2 = -f_bar / lambda_minor;
if a2 <= 0.0 || b2 <= 0.0 {
return None;
}
let a = a2.sqrt();
let b = b2.sqrt();
let vx = vec_major[0];
let vy = vec_major[1];
let mut phi = vy.atan2(vx);
if phi < 0.0 {
phi += std::f64::consts::PI;
}
Some(Ellipse {
center: Point::new(xc, yc),
semi_major: a,
semi_minor: b,
rotation: phi,
})
}
pub fn center(&self) -> Point {
self.center
}
pub fn semi_major(&self) -> f64 {
self.semi_major
}
pub fn semi_minor(&self) -> f64 {
self.semi_minor
}
pub fn rotation(&self) -> f64 {
self.rotation
}
pub fn sector_area(&self, theta: f64) -> f64 {
let a = self.semi_major;
let b = self.semi_minor;
let num = (b - a) * (2.0 * theta).sin();
let den = a + b + (b - a) * (2.0 * theta).cos();
0.5 * a * b * (theta - num.atan2(den))
}
pub fn sector_area_between(&self, theta0: f64, theta1: f64) -> f64 {
let t0 = theta0;
let mut t1 = theta1;
if t1 < t0 {
t1 += 2.0 * PI;
}
self.sector_area(t1) - self.sector_area(t0)
}
pub fn ellipse_segment(&self, p0: Point, p1: Point) -> f64 {
let p0 = p0.to_ellipse_frame(self);
let p1 = p1.to_ellipse_frame(self);
let theta0 = p0.angle_from_origin();
let mut theta1 = p1.angle_from_origin();
if theta1 < theta0 {
theta1 += 2.0 * PI;
}
let triangle = 0.5 * (p1.x() * p0.y() - p0.x() * p1.y()).abs();
if (theta1 - theta0) <= PI {
self.sector_area(theta1) - self.sector_area(theta0) - triangle
} else {
self.area() - (self.sector_area(theta0 + 2.0 * PI) - self.sector_area(theta1))
+ triangle
}
}
fn compute_lens_area(&self, other: &Self, p0: &Point, p1: &Point) -> f64 {
let seg1_10 = self.ellipse_segment(*p1, *p0);
let seg2_10 = other.ellipse_segment(*p1, *p0);
let triangle_10 = 0.5 * ((p0.x() + p1.x()) * (p0.y() - p1.y()));
let arc1 = triangle_10 + seg1_10.min(seg2_10);
let seg1_01 = self.ellipse_segment(*p0, *p1);
let seg2_01 = other.ellipse_segment(*p0, *p1);
let triangle_01 = 0.5 * ((p1.x() + p0.x()) * (p1.y() - p0.y()));
let arc2 = triangle_01 + seg1_01.min(seg2_01);
arc1 + arc2
}
fn compute_multi_point_intersection_area(&self, other: &Self, points: &[Point]) -> f64 {
if points.len() < 3 {
return 0.0;
}
let cx = points.iter().map(|p| p.x()).sum::<f64>() / points.len() as f64;
let cy = points.iter().map(|p| p.y()).sum::<f64>() / points.len() as f64;
let mut sorted_points: Vec<Point> = points.to_vec();
sorted_points.sort_by(|a, b| {
let angle_a = (a.x() - cx).atan2(a.y() - cy);
let angle_b = (b.x() - cx).atan2(b.y() - cy);
angle_a.partial_cmp(&angle_b).unwrap()
});
let mut total_area = 0.0;
let n = sorted_points.len();
for k in 0..n {
let l = if k == 0 { n - 1 } else { k - 1 };
let pi = sorted_points[k]; let pj = sorted_points[l];
let triangle = 0.5 * ((pj.x() + pi.x()) * (pj.y() - pi.y()));
let seg1 = self.ellipse_segment(pi, pj);
let seg2 = other.ellipse_segment(pi, pj);
total_area += triangle + seg1.min(seg2);
}
total_area
}
}
impl Area for Ellipse {
fn area(&self) -> f64 {
PI * self.semi_major * self.semi_minor
}
}
impl Perimeter for Ellipse {
fn perimeter(&self) -> f64 {
let a = self.semi_major;
let b = self.semi_minor;
let h = ((a - b).powi(2)) / ((a + b).powi(2));
PI * (a + b) * (1.0 + (3.0 * h) / (10.0 + (4.0 - 3.0 * h).sqrt()))
}
}
impl Centroid for Ellipse {
fn centroid(&self) -> Point {
self.center
}
}
impl BoundingBox for Ellipse {
fn bounding_box(&self) -> Rectangle {
let cos_theta = self.rotation.cos();
let sin_theta = self.rotation.sin();
let width = 2.0
* ((self.semi_major * cos_theta).powi(2) + (self.semi_minor * sin_theta).powi(2))
.sqrt();
let height = 2.0
* ((self.semi_major * sin_theta).powi(2) + (self.semi_minor * cos_theta).powi(2))
.sqrt();
Rectangle::new(self.center, width, height)
}
}
impl Closed for Ellipse {
fn contains(&self, other: &Self) -> bool {
if !self.contains_point(&other.center) {
return false;
}
if other.area() > self.area() {
return false;
}
let c1 = Conic::from_ellipse(*self);
let c2 = Conic::from_ellipse(*other);
let intersection_points = c1.intersect_conic(&c2);
if intersection_points.is_empty() {
return true;
}
let phi = other.rotation;
let cos_phi = phi.cos();
let sin_phi = phi.sin();
let major_offset_x = other.semi_major * cos_phi;
let major_offset_y = other.semi_major * sin_phi;
let p1 = Point::new(
other.center.x() + major_offset_x,
other.center.y() + major_offset_y,
);
let p2 = Point::new(
other.center.x() - major_offset_x,
other.center.y() - major_offset_y,
);
let minor_offset_x = other.semi_minor * -sin_phi;
let minor_offset_y = other.semi_minor * cos_phi;
let p3 = Point::new(
other.center.x() + minor_offset_x,
other.center.y() + minor_offset_y,
);
let p4 = Point::new(
other.center.x() - minor_offset_x,
other.center.y() - minor_offset_y,
);
self.contains_point(&p1)
&& self.contains_point(&p2)
&& self.contains_point(&p3)
&& self.contains_point(&p4)
}
fn contains_point(&self, point: &Point) -> bool {
let dx = point.x() - self.center.x();
let dy = point.y() - self.center.y();
let cos_phi = self.rotation.cos();
let sin_phi = self.rotation.sin();
let x_local = dx * cos_phi + dy * sin_phi;
let y_local = -dx * sin_phi + dy * cos_phi;
(x_local * x_local) / (self.semi_major * self.semi_major)
+ (y_local * y_local) / (self.semi_minor * self.semi_minor)
<= 1.0
}
fn intersects(&self, other: &Self) -> bool {
let center_distance = self.center.distance(&other.center);
let max_reach_self = self.semi_major.max(self.semi_minor);
let max_reach_other = other.semi_major.max(other.semi_minor);
if center_distance > max_reach_self + max_reach_other {
return false;
}
if self.contains(other) || other.contains(self) {
return false;
}
let c1 = Conic::from_ellipse(*self);
let c2 = Conic::from_ellipse(*other);
let intersection_points = c1.intersect_conic(&c2);
!intersection_points.is_empty()
}
fn intersection_area(&self, other: &Self) -> f64 {
if self.contains(other) {
return other.area();
}
if other.contains(self) {
return self.area();
}
let points = self.intersection_points(other);
let n = points.len();
match n {
0 => 0.0, 1 => {
0.0
}
2 => {
self.compute_lens_area(other, &points[0], &points[1])
}
_ => {
self.compute_multi_point_intersection_area(other, &points)
}
}
}
fn intersection_points(&self, other: &Self) -> Vec<Point> {
let c1 = Conic::from_ellipse(*self);
let c2 = Conic::from_ellipse(*other);
let homogeneous_points = c1.intersect_conic(&c2);
homogeneous_points
.into_iter()
.map(Point::from_homogeneous)
.collect()
}
}
#[derive(Debug, Clone, Copy)]
pub(crate) struct EllipseArc {
pub ellipse: usize,
pub t_start: f64,
pub delta_t: f64,
}
fn ellipse_point_at(e: &Ellipse, t: f64) -> Point {
let cos_t = t.cos();
let sin_t = t.sin();
let u = e.semi_major * cos_t;
let v = e.semi_minor * sin_t;
let cos_phi = e.rotation.cos();
let sin_phi = e.rotation.sin();
Point::new(
e.center.x() + u * cos_phi - v * sin_phi,
e.center.y() + u * sin_phi + v * cos_phi,
)
}
fn ellipse_param_for_point(e: &Ellipse, p: &Point) -> f64 {
let local = p.to_ellipse_frame(e);
(local.y() * e.semi_major).atan2(local.x() * e.semi_minor)
}
fn ellipse_implicit_value(p: &Point, e: &Ellipse) -> f64 {
let local = p.to_ellipse_frame(e);
(local.x() / e.semi_major).powi(2) + (local.y() / e.semi_minor).powi(2)
}
fn arc_midpoint_owned_by_j(j: usize, mid: &Point, indices: &[usize], ellipses: &[Ellipse]) -> bool {
for &l in indices {
if l == j {
continue;
}
let el = &ellipses[l];
let scale = el.semi_major.max(el.semi_minor);
let eps = (2.0 * BOUNDARY_COINCIDENCE_GEOM_TOL / scale).max(IMPLICIT_VALUE_FP_TOL);
let v = ellipse_implicit_value(mid, el);
if v > 1.0 + eps {
return false;
}
if (v - 1.0).abs() <= eps && l < j {
return false;
}
}
true
}
const BOUNDARY_COINCIDENCE_GEOM_TOL: f64 = 1e-7;
const IMPLICIT_VALUE_FP_TOL: f64 = 1e-12;
pub(crate) fn region_boundary_arcs_ellipse(
mask: RegionMask,
ellipses: &[Ellipse],
intersections: &[crate::geometry::diagram::IntersectionPoint],
n_sets: usize,
) -> Vec<EllipseArc> {
use crate::geometry::diagram::{adopters_to_mask, mask_to_indices};
let count = mask.count_ones();
if count == 0 {
return Vec::new();
}
if count == 1 {
let idx = mask.trailing_zeros() as usize;
return vec![EllipseArc {
ellipse: idx,
t_start: 0.0,
delta_t: 2.0 * PI,
}];
}
let indices = mask_to_indices(mask, n_sets);
let mut arcs: Vec<EllipseArc> = Vec::new();
let two_pi = 2.0 * PI;
let arc_eps = 1e-9;
for &j in &indices {
let ej = &ellipses[j];
let mut j_ts: Vec<f64> = intersections
.iter()
.filter(|ip| {
let (p1, p2) = ip.parents();
if p1 != j && p2 != j {
return false;
}
let am = adopters_to_mask(ip.adopters());
(mask & am) == mask
})
.map(|ip| ellipse_param_for_point(ej, ip.point()))
.collect();
if j_ts.is_empty() {
let probe = ellipse_point_at(ej, 0.0);
if arc_midpoint_owned_by_j(j, &probe, &indices, ellipses) {
arcs.push(EllipseArc {
ellipse: j,
t_start: 0.0,
delta_t: two_pi,
});
}
continue;
}
j_ts.sort_by(|a, b| a.partial_cmp(b).unwrap_or(std::cmp::Ordering::Equal));
let m = j_ts.len();
for k in 0..m {
let t_a = j_ts[k];
let t_b = j_ts[(k + 1) % m];
let mut delta = t_b - t_a;
if delta <= 0.0 {
delta += two_pi;
}
if delta < arc_eps || delta > two_pi - arc_eps {
continue;
}
let t_mid = t_a + delta * 0.5;
let mid = ellipse_point_at(ej, t_mid);
if !arc_midpoint_owned_by_j(j, &mid, &indices, ellipses) {
continue;
}
arcs.push(EllipseArc {
ellipse: j,
t_start: t_a,
delta_t: delta,
});
}
}
arcs
}
pub(crate) fn area_from_boundary_arcs_ellipse(arcs: &[EllipseArc], ellipses: &[Ellipse]) -> f64 {
let mut total = 0.0;
for arc in arcs {
let e = &ellipses[arc.ellipse];
let xc = e.center.x();
let yc = e.center.y();
let a = e.semi_major;
let b = e.semi_minor;
let cphi = e.rotation.cos();
let sphi = e.rotation.sin();
let t_a = arc.t_start;
let t_b = t_a + arc.delta_t;
let s_b_minus_s_a = t_b.sin() - t_a.sin();
let c_b_minus_c_a = t_b.cos() - t_a.cos();
let int_x = b * cphi * s_b_minus_s_a + a * sphi * c_b_minus_c_a;
let int_y = -a * cphi * c_b_minus_c_a + b * sphi * s_b_minus_s_a;
total += xc * int_x + yc * int_y + a * b * arc.delta_t;
}
0.5 * total
}
pub(crate) fn accumulate_region_overlap_gradient_ellipse(
arcs: &[EllipseArc],
ellipses: &[Ellipse],
grad: &mut [f64],
) {
for arc in arcs {
let e = &ellipses[arc.ellipse];
let a = e.semi_major;
let b = e.semi_minor;
let cphi = e.rotation.cos();
let sphi = e.rotation.sin();
let t_a = arc.t_start;
let t_b = t_a + arc.delta_t;
let s_a = t_a.sin();
let s_b = t_b.sin();
let c_a = t_a.cos();
let c_b = t_b.cos();
let s2_a = (2.0 * t_a).sin();
let s2_b = (2.0 * t_b).sin();
let c2_a = (2.0 * t_a).cos();
let c2_b = (2.0 * t_b).cos();
let off = arc.ellipse * 5;
grad[off] += b * cphi * (s_b - s_a) + a * sphi * (c_b - c_a);
grad[off + 1] += -a * cphi * (c_b - c_a) + b * sphi * (s_b - s_a);
grad[off + 2] += a * (0.5 * b * arc.delta_t + 0.25 * b * (s2_b - s2_a));
grad[off + 3] += b * (0.5 * a * arc.delta_t - 0.25 * a * (s2_b - s2_a));
grad[off + 4] += 0.25 * (a * a - b * b) * (c2_a - c2_b);
}
}
pub(crate) fn collect_intersections_ellipse(
ellipses: &[Ellipse],
) -> Vec<crate::geometry::diagram::IntersectionPoint> {
use crate::geometry::diagram::IntersectionPoint;
let n_sets = ellipses.len();
let mut intersections: Vec<IntersectionPoint> = Vec::new();
for i in 0..n_sets {
for j in (i + 1)..n_sets {
let pts = ellipses[i].intersection_points(&ellipses[j]);
for p in pts {
let adopters: Vec<usize> = (0..n_sets)
.filter(|&k| {
let local = p.to_ellipse_frame(&ellipses[k]);
let val = (local.x() / ellipses[k].semi_major).powi(2)
+ (local.y() / ellipses[k].semi_minor).powi(2);
val <= 1.0 + 1e-9
})
.collect();
intersections.push(IntersectionPoint::new(p, (i, j), adopters));
}
}
}
intersections
}
impl DiagramShape for Ellipse {
fn compute_exclusive_regions(shapes: &[Self]) -> std::collections::HashMap<RegionMask, f64> {
use crate::geometry::diagram::{discover_regions, to_exclusive_areas};
use std::collections::HashMap;
let n_sets = shapes.len();
let intersections = collect_intersections_ellipse(shapes);
let regions = discover_regions(shapes, &intersections, n_sets);
let mut overlapping_areas: HashMap<RegionMask, f64> = HashMap::new();
for &mask in ®ions {
let arcs = region_boundary_arcs_ellipse(mask, shapes, &intersections, n_sets);
overlapping_areas.insert(mask, area_from_boundary_arcs_ellipse(&arcs, shapes));
}
to_exclusive_areas(&overlapping_areas)
}
fn compute_exclusive_regions_with_gradient(
shapes: &[Self],
) -> Option<crate::geometry::traits::ExclusiveRegionsAndGradient> {
use crate::geometry::diagram::{discover_regions, to_exclusive_areas_and_gradients};
use std::collections::HashMap;
let n_sets = shapes.len();
let n_params = n_sets * 5;
let intersections = collect_intersections_ellipse(shapes);
let regions = discover_regions(shapes, &intersections, n_sets);
let mut overlapping_areas: HashMap<RegionMask, f64> = HashMap::new();
let mut overlapping_grads: HashMap<RegionMask, Vec<f64>> = HashMap::new();
for &mask in ®ions {
let arcs = region_boundary_arcs_ellipse(mask, shapes, &intersections, n_sets);
overlapping_areas.insert(mask, area_from_boundary_arcs_ellipse(&arcs, shapes));
let mut grad = vec![0.0; n_params];
accumulate_region_overlap_gradient_ellipse(&arcs, shapes, &mut grad);
overlapping_grads.insert(mask, grad);
}
Some(to_exclusive_areas_and_gradients(
&overlapping_areas,
&overlapping_grads,
n_params,
))
}
fn optimizer_params_from_circle(x: f64, y: f64, radius: f64) -> Vec<f64> {
let log_radius = radius.ln();
vec![x, y, log_radius, log_radius, 0.0]
}
fn mds_target_distance(
area_i: f64,
area_j: f64,
target_overlap: f64,
) -> Result<f64, crate::error::DiagramError> {
Circle::mds_target_distance(area_i, area_j, target_overlap)
}
fn n_params() -> usize {
5 }
fn from_params(params: &[f64]) -> Self {
assert_eq!(
params.len(),
5,
"Ellipse requires 5 parameters: x, y, semi_major, semi_minor, rotation"
);
Ellipse::new(
Point::new(params[0], params[1]),
params[2],
params[3],
params[4],
)
}
fn to_params(&self) -> Vec<f64> {
vec![
self.center.x(),
self.center.y(),
self.semi_major,
self.semi_minor,
self.rotation,
]
}
fn from_optimizer_params(params: &[f64]) -> Self {
assert_eq!(
params.len(),
5,
"Ellipse requires 5 parameters: x, y, ln(semi_major), ln(semi_minor), rotation"
);
Ellipse::new(
Point::new(params[0], params[1]),
params[2].exp().max(f64::MIN_POSITIVE),
params[3].exp().max(f64::MIN_POSITIVE),
params[4],
)
}
fn to_optimizer_params(&self) -> Vec<f64> {
vec![
self.center.x(),
self.center.y(),
self.semi_major.ln(),
self.semi_minor.ln(),
self.rotation,
]
}
fn canonical_venn_layout(n: usize) -> Option<Vec<Self>> {
let params: &[CanonicalVennParams] = match n {
1 => &VENN_N1,
2 => &VENN_N2,
3 => &VENN_N3,
4 => &VENN_N4,
5 => &VENN_N5,
_ => return None,
};
Some(
params
.iter()
.map(|&(h, k, a, b, phi)| Ellipse::new(Point::new(h, k), a, b, phi))
.collect(),
)
}
}
type CanonicalVennParams = (f64, f64, f64, f64, f64);
const VENN_N1: [CanonicalVennParams; 1] = [(0.0, 0.0, 1.0, 1.0, 0.0)];
const VENN_N2: [CanonicalVennParams; 2] = [(-0.5, 0.0, 1.0, 1.0, 1.0), (0.5, 0.0, 1.0, 1.0, 1.0)];
const VENN_N3: [CanonicalVennParams; 3] = [
(-0.42, -0.36, 1.05, 1.05, 3.76),
(0.42, -0.36, 1.05, 1.05, 3.76),
(0.00, 0.36, 1.05, 1.05, 3.76),
];
const VENN_N4: [CanonicalVennParams; 4] = [
(-0.8, 0.0, 1.2, 2.0, PI / 4.0),
(0.8, 0.0, 1.2, 2.0, -PI / 4.0),
(0.0, 1.0, 1.2, 2.0, PI / 4.0),
(0.0, 1.0, 1.2, 2.0, -PI / 4.0),
];
const VENN_N5: [CanonicalVennParams; 5] = [
(0.176, 0.096, 1.0, 0.6, 0.000),
(-0.037, 0.197, 1.0, 0.6, 1.257),
(-0.198, 0.026, 1.0, 0.6, 2.513),
(-0.086, -0.181, 1.0, 0.6, 3.770),
(0.145, -0.137, 1.0, 0.6, 5.027),
];
impl Polygonize for Ellipse {
fn polygonize(&self, n_vertices: usize) -> Polygon {
use std::f64::consts::PI;
let n_vertices = n_vertices.max(3);
let mut vertices = Vec::with_capacity(n_vertices);
let cos_rot = self.rotation.cos();
let sin_rot = self.rotation.sin();
for i in 0..n_vertices {
let angle = 2.0 * PI * (i as f64) / (n_vertices as f64);
let local_x = self.semi_major * angle.cos();
let local_y = self.semi_minor * angle.sin();
let x = self.center.x() + local_x * cos_rot - local_y * sin_rot;
let y = self.center.y() + local_x * sin_rot + local_y * cos_rot;
vertices.push(Point::new(x, y));
}
Polygon::new(vertices)
}
}
#[cfg(test)]
mod tests {
use super::*;
use crate::geometry::shapes::Circle;
fn approx_eq(a: f64, b: f64) -> bool {
(a - b).abs() < 1e-10
}
#[test]
fn test_sector_area_between_quarter() {
let ellipse = Ellipse::new(Point::new(0.0, 0.0), 5.0, 3.0, 0.0);
let sector_between = ellipse.sector_area_between(0.0, PI / 2.0);
let expected = ellipse.sector_area(PI / 2.0) - ellipse.sector_area(0.0);
assert!(approx_eq(sector_between, expected));
assert!(approx_eq(sector_between, ellipse.area() / 4.0));
}
#[test]
fn test_sector_area_between_with_ccw_adjustment() {
let ellipse = Ellipse::new(Point::new(0.0, 0.0), 4.0, 2.0, 0.0);
let theta0 = 3.0 * PI / 4.0;
let theta1 = PI / 4.0;
let sector_between = ellipse.sector_area_between(theta0, theta1);
let expected = ellipse.sector_area(theta1 + 2.0 * PI) - ellipse.sector_area(theta0);
assert!(
approx_eq(sector_between, expected),
"Expected {}, got {}",
expected,
sector_between
);
assert!(sector_between > 0.0);
assert!(sector_between > ellipse.area() / 2.0);
}
#[test]
fn test_sector_area_between_same_angle() {
let ellipse = Ellipse::new(Point::new(1.0, 2.0), 3.0, 2.5, 0.3);
let angle = PI / 3.0;
let sector = ellipse.sector_area_between(angle, angle);
assert!(approx_eq(sector, 0.0));
}
#[test]
fn test_sector_area_between_full_rotation() {
let ellipse = Ellipse::new(Point::new(0.0, 0.0), 5.0, 3.0, 0.0);
let sector = ellipse.sector_area_between(0.0, 2.0 * PI);
assert!(approx_eq(sector, ellipse.area()));
}
#[test]
fn test_sector_full_ellipse() {
let ellipse = Ellipse::new(Point::new(0.0, 0.0), 3.0, 2.0, 0.0);
let full_angle = 2.0 * PI;
let sector = ellipse.sector_area(full_angle);
assert!(approx_eq(sector, ellipse.area()));
}
#[test]
fn test_sector_half_ellipse() {
let ellipse = Ellipse::new(Point::new(0.0, 0.0), 4.0, 3.0, 0.0);
let half_angle = PI;
let sector = ellipse.sector_area(half_angle);
assert!(approx_eq(sector, ellipse.area() / 2.0));
}
#[test]
fn test_sector_quarter_ellipse() {
let ellipse = Ellipse::new(Point::new(0.0, 0.0), 5.0, 3.0, 0.0);
let quarter_angle = PI / 2.0;
let sector = ellipse.sector_area(quarter_angle);
assert!(approx_eq(sector, ellipse.area() / 4.0));
}
#[test]
fn test_sector_zero_angle() {
let ellipse = Ellipse::new(Point::new(1.0, 1.0), 3.0, 2.0, 0.5);
let sector = ellipse.sector_area(0.0);
assert!(approx_eq(sector, 0.0));
}
#[test]
fn test_sector_circular_ellipse_matches_circle() {
let radius = 3.0;
let ellipse = Ellipse::new(Point::new(0.0, 0.0), radius, radius, 0.0);
let circle = Circle::new(Point::new(0.0, 0.0), radius);
let angle = PI / 3.0;
let ellipse_sector = ellipse.sector_area(angle);
let circle_sector = circle.sector_area(angle);
assert!(approx_eq(ellipse_sector, circle_sector));
}
#[test]
fn test_sector_circular_ellipse_multiple_angles() {
let radius = 2.0;
let ellipse = Ellipse::new(Point::new(0.0, 0.0), radius, radius, 0.0);
let circle = Circle::new(Point::new(0.0, 0.0), radius);
let test_angles = vec![
0.0,
PI / 6.0,
PI / 4.0,
PI / 2.0,
PI,
3.0 * PI / 2.0,
2.0 * PI,
];
for angle in test_angles {
let ellipse_sector = ellipse.sector_area(angle);
let circle_sector = circle.sector_area(angle);
assert!(
approx_eq(ellipse_sector, circle_sector),
"Mismatch at angle {}: ellipse={}, circle={}",
angle,
ellipse_sector,
circle_sector
);
}
}
#[test]
fn test_ellipse_segment_semicircle() {
let radius = 2.0;
let ellipse = Ellipse::new(Point::new(0.0, 0.0), radius, radius, 0.0);
let p0 = Point::new(-radius, 0.0);
let p1 = Point::new(radius, 0.0);
let segment = ellipse.ellipse_segment(p0, p1);
assert!(approx_eq(segment, ellipse.area() / 2.0));
}
#[test]
fn test_ellipse_segment_circular_matches_circle_segment() {
let radius = 3.0;
let ellipse = Ellipse::new(Point::new(0.0, 0.0), radius, radius, 0.0);
let circle = Circle::new(Point::new(0.0, 0.0), radius);
let angle = PI / 4.0;
let p0 = Point::new(radius, 0.0);
let p1 = Point::new(radius * angle.cos(), radius * angle.sin());
let ellipse_segment = ellipse.ellipse_segment(p0, p1);
let circle_segment = circle.segment_area_from_angle(angle);
assert!(
approx_eq(ellipse_segment, circle_segment),
"Circular ellipse segment should match circle segment: ellipse={}, circle={}",
ellipse_segment,
circle_segment
);
}
#[test]
fn test_ellipse_segment_small_angle() {
let ellipse = Ellipse::new(Point::new(0.0, 0.0), 3.0, 2.0, 0.0);
let angle = PI / 12.0; let p0 = Point::new(ellipse.semi_major(), 0.0);
let p1 = Point::new(
ellipse.semi_major() * angle.cos(),
ellipse.semi_minor() * angle.sin(),
);
let segment = ellipse.ellipse_segment(p0, p1);
assert!(segment > 0.0);
let sector = ellipse.sector_area(angle);
assert!(segment < sector);
}
#[test]
fn test_ellipse_segment_zero_angle() {
let ellipse = Ellipse::new(Point::new(0.0, 0.0), 4.0, 2.5, 0.0);
let p = Point::new(ellipse.semi_major(), 0.0);
let segment = ellipse.ellipse_segment(p, p);
assert!(approx_eq(segment, 0.0));
}
#[test]
fn test_ellipse_segment_with_known_values() {
let ellipse = Ellipse::new(Point::new(0.0, 0.0), 5.0, 3.0, 0.0);
let p0 = Point::new(ellipse.semi_major(), 0.0);
let p1 = Point::new(0.0, ellipse.semi_minor());
let segment_quarter = ellipse.ellipse_segment(p0, p1);
let expected_quarter = 4.280972450961723;
assert!(
approx_eq(segment_quarter, expected_quarter),
"Quarter segment: expected {}, got {}",
expected_quarter,
segment_quarter
);
let p2 = Point::new(-ellipse.semi_major(), 0.0);
let segment_half = ellipse.ellipse_segment(p0, p2);
let expected_half = 23.561944901923447;
assert!(
approx_eq(segment_half, expected_half),
"Half segment: expected {}, got {}",
expected_half,
segment_half
);
assert!(approx_eq(segment_half, ellipse.area() / 2.0));
}
#[test]
fn test_ellipse_segment_additional_angles() {
let ellipse = Ellipse::new(Point::new(0.0, 0.0), 5.0, 3.0, 0.0);
let angle = PI / 6.0;
let p0 = Point::new(ellipse.semi_major(), 0.0);
let p1 = Point::new(
ellipse.semi_major() * angle.cos(),
ellipse.semi_minor() * angle.sin(),
);
let segment_30deg = ellipse.ellipse_segment(p0, p1);
let expected_30deg = 0.17699081698724095;
assert!(
approx_eq(segment_30deg, expected_30deg),
"30° segment: expected {}, got {}",
expected_30deg,
segment_30deg
);
let angle = PI / 3.0;
let p2 = Point::new(
ellipse.semi_major() * angle.cos(),
ellipse.semi_minor() * angle.sin(),
);
let segment_60deg = ellipse.ellipse_segment(p0, p2);
let expected_60deg = 1.3587911055911919;
assert!(
approx_eq(segment_60deg, expected_60deg),
"60° segment: expected {}, got {}",
expected_60deg,
segment_60deg
);
assert!(segment_30deg < segment_60deg);
assert!(segment_60deg < ellipse.area() / 4.0);
}
#[test]
fn test_ellipse_segment_270_degrees() {
let ellipse = Ellipse::new(Point::new(0.0, 0.0), 5.0, 3.0, 0.0);
let angle = 3.0 * PI / 2.0;
let p0 = Point::new(ellipse.semi_major(), 0.0);
let p1 = Point::new(
ellipse.semi_major() * angle.cos(),
ellipse.semi_minor() * angle.sin(),
);
let segment_270deg = ellipse.ellipse_segment(p0, p1);
let expected_270deg = 42.84291735288517;
assert!(
approx_eq(segment_270deg, expected_270deg),
"270° segment: expected {}, got {}",
expected_270deg,
segment_270deg
);
assert!(segment_270deg > ellipse.area() / 2.0);
let complementary = ellipse.area() - segment_270deg;
assert!(complementary > 0.0);
assert!(complementary < ellipse.area() / 4.0);
}
#[test]
fn test_contains_point_center_and_boundary() {
let ellipse = Ellipse::new(Point::new(2.0, 3.0), 5.0, 3.0, 0.0);
assert!(ellipse.contains_point(&Point::new(2.0, 3.0)));
assert!(ellipse.contains_point(&Point::new(7.0, 3.0))); assert!(ellipse.contains_point(&Point::new(-3.0, 3.0)));
assert!(ellipse.contains_point(&Point::new(2.0, 6.0))); assert!(ellipse.contains_point(&Point::new(2.0, 0.0)));
assert!(!ellipse.contains_point(&Point::new(10.0, 10.0)));
}
#[test]
fn test_contains_point_rotated_ellipse() {
let ellipse = Ellipse::new(Point::new(0.0, 0.0), 4.0, 2.0, PI / 4.0);
assert!(ellipse.contains_point(&Point::new(0.0, 0.0)));
let dist = 4.0 / (2.0_f64).sqrt(); assert!(ellipse.contains_point(&Point::new(dist, dist)));
assert!(ellipse.contains_point(&Point::new(-dist, -dist)));
assert!(!ellipse.contains_point(&Point::new(5.0, 0.0)));
assert!(!ellipse.contains_point(&Point::new(0.0, 5.0)));
assert!(ellipse.contains_point(&Point::new(0.5, 0.5)));
}
#[test]
fn test_contains_ellipse_concentric() {
let outer = Ellipse::new(Point::new(0.0, 0.0), 5.0, 3.0, 0.0);
let inner = Ellipse::new(Point::new(0.0, 0.0), 3.0, 2.0, 0.0);
assert!(outer.contains(&inner));
assert!(!inner.contains(&outer));
}
#[test]
fn test_contains_ellipse_offset_centers() {
let outer = Ellipse::new(Point::new(0.0, 0.0), 10.0, 8.0, 0.0);
let inner = Ellipse::new(Point::new(2.0, 1.0), 2.0, 1.5, 0.0);
assert!(outer.contains(&inner));
let outside = Ellipse::new(Point::new(9.0, 0.0), 2.0, 1.5, 0.0);
assert!(!outer.contains(&outside));
}
#[test]
fn test_contains_ellipse_rotated() {
let outer = Ellipse::new(Point::new(0.0, 0.0), 6.0, 4.0, PI / 6.0);
let inner = Ellipse::new(Point::new(0.0, 0.0), 2.0, 1.5, PI / 4.0);
assert!(outer.contains(&inner));
assert!(!inner.contains(&outer));
}
#[test]
fn test_contains_ellipse_intersecting() {
let e1 = Ellipse::new(Point::new(0.0, 0.0), 5.0, 3.0, 0.0);
let e2 = Ellipse::new(Point::new(3.0, 0.0), 4.0, 2.5, 0.0);
assert!(!e1.contains(&e2));
assert!(!e2.contains(&e1));
}
#[test]
fn test_intersects_overlapping() {
let e1 = Ellipse::new(Point::new(0.0, 0.0), 5.0, 3.0, 0.0);
let e2 = Ellipse::new(Point::new(3.0, 0.0), 4.0, 2.5, 0.0);
assert!(e1.intersects(&e2));
assert!(e2.intersects(&e1));
}
#[test]
fn test_intersects_separate() {
let e1 = Ellipse::new(Point::new(0.0, 0.0), 2.0, 1.5, 0.0);
let e2 = Ellipse::new(Point::new(10.0, 0.0), 3.0, 2.0, 0.0);
assert!(!e1.intersects(&e2));
assert!(!e2.intersects(&e1));
}
#[test]
fn test_intersects_contained() {
let outer = Ellipse::new(Point::new(0.0, 0.0), 5.0, 3.0, 0.0);
let inner = Ellipse::new(Point::new(0.0, 0.0), 3.0, 2.0, 0.0);
assert!(!outer.intersects(&inner));
assert!(!inner.intersects(&outer));
}
#[test]
fn test_intersects_touching() {
let e1 = Ellipse::new(Point::new(0.0, 0.0), 3.0, 2.0, 0.0);
let e2 = Ellipse::new(Point::new(6.0, 0.0), 3.0, 2.0, 0.0);
assert!(e1.intersects(&e2));
}
#[test]
fn test_intersects_when_semi_minor_greater_than_semi_major() {
let e0 = Ellipse::new(
Point::new(36.86134450095962, 78.5649458016757),
(2.9862175109548414f64).exp(),
(4.125140737010668f64).exp(),
0.01270551956797799,
);
let e1 = Ellipse::new(
Point::new(86.13182443212173, 20.46219955112993),
(4.000229123217676f64).exp(),
(3.0943704188347185f64).exp(),
-0.023783971647556423,
);
assert!(
e0.intersection_area(&e1) > 0.0,
"sanity: ellipses do overlap by area",
);
assert!(
e0.intersects(&e1),
"intersects must agree with intersection_area on overlap",
);
}
#[test]
fn test_intersection_points_overlapping() {
let e1 = Ellipse::new(Point::new(0.0, 0.0), 5.0, 3.0, 0.0);
let e2 = Ellipse::new(Point::new(3.0, 0.0), 4.0, 2.5, 0.0);
let points = e1.intersection_points(&e2);
assert!(!points.is_empty());
assert!(points.len() <= 4);
}
#[test]
fn test_intersection_points_separate() {
let e1 = Ellipse::new(Point::new(0.0, 0.0), 2.0, 1.5, 0.0);
let e2 = Ellipse::new(Point::new(10.0, 0.0), 3.0, 2.0, 0.0);
let points = e1.intersection_points(&e2);
assert!(points.len() <= 4);
}
#[test]
fn test_intersection_points_contained() {
let outer = Ellipse::new(Point::new(0.0, 0.0), 5.0, 3.0, 0.0);
let inner = Ellipse::new(Point::new(0.0, 0.0), 3.0, 2.0, 0.0);
let points = outer.intersection_points(&inner);
assert!(points.len() <= 4);
}
#[test]
fn test_intersection_points_tangent() {
let e1 = Ellipse::new(Point::new(0.0, 0.0), 2.0, 2.0, 0.0); let e2 = Ellipse::new(Point::new(4.0, 0.0), 2.0, 2.0, 0.0);
let points = e1.intersection_points(&e2);
assert!(!points.is_empty());
if !points.is_empty() {
let first_point = &points[0];
assert!((first_point.x() - 2.0).abs() < 1e-6);
assert!(first_point.y().abs() < 1e-6);
}
}
#[test]
fn test_intersection_area_no_overlap() {
let e1 = Ellipse::new(Point::new(0.0, 0.0), 2.0, 1.5, 0.0);
let e2 = Ellipse::new(Point::new(10.0, 0.0), 2.0, 1.5, 0.0);
let area = e1.intersection_area(&e2);
assert!(approx_eq(area, 0.0), "Expected 0, got {}", area);
}
#[test]
fn test_intersection_area_one_contains_other() {
let e1 = Ellipse::new(Point::new(0.0, 0.0), 5.0, 4.0, 0.0);
let e2 = Ellipse::new(Point::new(0.0, 0.0), 2.0, 1.5, 0.0);
let area = e1.intersection_area(&e2);
let expected = e2.area();
assert!(
(area - expected).abs() < 1e-8,
"Expected {}, got {}",
expected,
area
);
let area_reverse = e2.intersection_area(&e1);
assert!(
(area_reverse - expected).abs() < 1e-8,
"Expected {}, got {}",
expected,
area_reverse
);
}
#[test]
fn test_intersection_area_identical_ellipses() {
let e1 = Ellipse::new(Point::new(1.0, 2.0), 3.0, 2.0, PI / 4.0);
let e2 = Ellipse::new(Point::new(1.0, 2.0), 3.0, 2.0, PI / 4.0);
let area = e1.intersection_area(&e2);
let expected = e1.area();
assert!(
(area - expected).abs() < 1e-8,
"Expected {}, got {}",
expected,
area
);
}
#[test]
fn test_intersection_area_circles_partial_overlap() {
let radius = 3.0;
let e1 = Ellipse::new(Point::new(0.0, 0.0), radius, radius, 0.0);
let e2 = Ellipse::new(Point::new(4.0, 0.0), radius, radius, 0.0);
let area = e1.intersection_area(&e2);
let c1 = Circle::new(Point::new(0.0, 0.0), radius);
let c2 = Circle::new(Point::new(4.0, 0.0), radius);
let expected = c1.intersection_area(&c2);
assert!(
(area - expected).abs() < 1e-6,
"Expected {}, got {}",
expected,
area
);
}
#[test]
fn test_intersection_area_tangent_ellipses() {
let e1 = Ellipse::new(Point::new(0.0, 0.0), 2.0, 2.0, 0.0);
let e2 = Ellipse::new(Point::new(4.0, 0.0), 2.0, 2.0, 0.0);
let area = e1.intersection_area(&e2);
assert!(area < 1e-6, "Expected ~0, got {}", area);
}
#[test]
fn test_intersection_area_symmetric() {
let e1 = Ellipse::new(Point::new(0.0, 0.0), 4.0, 3.0, 0.0);
let e2 = Ellipse::new(Point::new(3.0, 2.0), 3.5, 2.5, PI / 6.0);
let area1 = e1.intersection_area(&e2);
let area2 = e2.intersection_area(&e1);
assert!(
(area1 - area2).abs() < 1e-8,
"Areas should be symmetric: {} vs {}",
area1,
area2
);
}
#[test]
fn test_intersection_area_bounds() {
let e1 = Ellipse::new(Point::new(0.0, 0.0), 5.0, 4.0, 0.0);
let e2 = Ellipse::new(Point::new(2.0, 1.0), 3.0, 2.5, PI / 8.0);
let area = e1.intersection_area(&e2);
let min_area = e1.area().min(e2.area());
assert!(
area >= 0.0,
"Intersection area should be non-negative: {}",
area
);
assert!(
area <= min_area + 1e-6,
"Intersection area {} should not exceed smaller ellipse area {}",
area,
min_area
);
}
#[test]
fn test_intersection_area_rotated_ellipses() {
let e1 = Ellipse::new(Point::new(0.0, 0.0), 4.0, 2.5, 0.0);
let e2 = Ellipse::new(Point::new(3.0, 0.0), 4.0, 2.5, PI / 2.0);
let area = e1.intersection_area(&e2);
assert!(area > 0.0, "Expected positive area, got {}", area);
assert!(
area < e1.area() && area < e2.area(),
"Area {} should be less than both ellipse areas",
area
);
}
#[test]
fn test_intersection_area_nearly_identical() {
let e1 = Ellipse::new(Point::new(0.0, 0.0), 3.0, 2.0, 0.0);
let e2 = Ellipse::new(Point::new(0.01, 0.01), 3.0, 2.0, 0.0);
let area = e1.intersection_area(&e2);
let expected = e1.area();
assert!(
(area - expected).abs() < 0.1,
"Expected ~{}, got {}",
expected,
area
);
}
fn monte_carlo_intersection_area(
e1: &Ellipse,
e2: &Ellipse,
n_samples: usize,
seed: u64,
) -> f64 {
use rand::rngs::StdRng;
use rand::{Rng, SeedableRng};
let mut rng = StdRng::seed_from_u64(seed);
let mut in_both = 0;
let bbox1 = e1.bounding_box();
let bbox2 = e2.bounding_box();
let (min1, max1) = bbox1.to_points();
let (min2, max2) = bbox2.to_points();
let x_min = min1.x().max(min2.x());
let x_max = max1.x().min(max2.x());
let y_min = min1.y().max(min2.y());
let y_max = max1.y().min(max2.y());
if x_min >= x_max || y_min >= y_max {
return 0.0; }
let bbox_area = (x_max - x_min) * (y_max - y_min);
for _ in 0..n_samples {
let x = x_min + (x_max - x_min) * rng.random::<f64>();
let y = y_min + (y_max - y_min) * rng.random::<f64>();
let p = Point::new(x, y);
if e1.contains_point(&p) && e2.contains_point(&p) {
in_both += 1;
}
}
(in_both as f64 / n_samples as f64) * bbox_area
}
#[test]
fn test_intersection_area_monte_carlo_circles() {
let radius = 3.0;
let e1 = Ellipse::new(Point::new(0.0, 0.0), radius, radius, 0.0);
let e2 = Ellipse::new(Point::new(4.0, 0.0), radius, radius, 0.0);
let exact = e1.intersection_area(&e2);
let mc_estimate = monte_carlo_intersection_area(&e1, &e2, 100_000, 42);
let error = (exact - mc_estimate).abs() / exact;
assert!(
error < 0.05,
"Monte Carlo and exact should agree within 5%: exact={}, mc={}, error={:.1}%",
exact,
mc_estimate,
error * 100.0
);
}
#[test]
fn test_intersection_area_monte_carlo_rotated_ellipses() {
let e1 = Ellipse::new(Point::new(0.0, 0.0), 4.0, 2.5, 0.0);
let e2 = Ellipse::new(Point::new(3.0, 0.0), 4.0, 2.5, PI / 2.0);
let exact = e1.intersection_area(&e2);
let mc_estimate = monte_carlo_intersection_area(&e1, &e2, 100_000, 42);
let error = (exact - mc_estimate).abs() / exact.max(mc_estimate);
assert!(
error < 0.05,
"Monte Carlo and exact should agree within 5%: exact={}, mc={}, error={:.1}%",
exact,
mc_estimate,
error * 100.0
);
}
#[test]
fn test_intersection_area_monte_carlo_general_ellipses() {
let e1 = Ellipse::new(Point::new(1.0, 0.5), 3.5, 2.0, PI / 6.0);
let e2 = Ellipse::new(Point::new(3.0, 1.5), 3.0, 2.5, PI / 4.0);
let exact = e1.intersection_area(&e2);
let mc_estimate = monte_carlo_intersection_area(&e1, &e2, 150_000, 42);
if exact > 0.5 && mc_estimate > 0.5 {
let error = (exact - mc_estimate).abs() / exact.max(mc_estimate);
assert!(
error < 0.06,
"Monte Carlo and exact should agree within 6%: exact={}, mc={}, error={:.1}%",
exact,
mc_estimate,
error * 100.0
);
} else {
assert!(
exact < 1.0 && mc_estimate < 1.0,
"Both should be small: exact={}, mc={}",
exact,
mc_estimate
);
}
}
#[test]
fn test_intersection_area_monte_carlo_contained() {
let e1 = Ellipse::new(Point::new(0.0, 0.0), 5.0, 4.0, 0.0);
let e2 = Ellipse::new(Point::new(0.5, 0.5), 2.0, 1.5, PI / 8.0);
let exact = e1.intersection_area(&e2);
let expected = e2.area();
assert!(e1.contains(&e2), "e1 should contain e2");
assert!(
(exact - expected).abs() < 1e-6,
"Intersection should equal smaller ellipse area: exact={}, expected={}",
exact,
expected
);
let mc_estimate = monte_carlo_intersection_area(&e1, &e2, 100_000, 42);
let error = (exact - mc_estimate).abs() / exact;
assert!(
error < 0.03,
"Monte Carlo should agree with exact (contained case): exact={}, mc={}, error={:.1}%",
exact,
mc_estimate,
error * 100.0
);
}
#[test]
fn test_from_radius_ratio_circle() {
let e = Ellipse::from_radius_ratio(Point::new(1.0, 2.0), 3.0, 1.0, 0.5);
assert!((e.semi_major() - 3.0).abs() < 1e-10);
assert!((e.semi_minor() - 3.0).abs() < 1e-10);
assert_eq!(e.center(), Point::new(1.0, 2.0));
assert_eq!(e.rotation(), 0.5);
}
#[test]
fn test_from_radius_ratio_elongated() {
let radius = 2.0;
let log_aspect = -std::f64::consts::LN_2; let e = Ellipse::from_radius_ratio(Point::new(0.0, 0.0), radius, log_aspect, 0.0);
let aspect = log_aspect.exp();
let expected_a = radius / aspect.sqrt();
let expected_b = radius * aspect.sqrt();
assert!((e.semi_major() - expected_a).abs() < 1e-10);
assert!((e.semi_minor() - expected_b).abs() < 1e-10);
let geom_mean = (e.semi_major() * e.semi_minor()).sqrt();
assert!((geom_mean - radius).abs() < 1e-10);
}
#[test]
fn test_from_radius_ratio_preserves_area() {
let radius = 3.0;
let expected_area = PI * radius * radius;
for log_aspect in [-1.6, -0.693, -0.223, 0.0] {
let e = Ellipse::from_radius_ratio(Point::new(0.0, 0.0), radius, log_aspect, 0.0);
assert!((e.area() - expected_area).abs() < 1e-9);
}
}
#[test]
fn test_from_radius_ratio_clamping() {
let e_low = Ellipse::from_radius_ratio(Point::new(0.0, 0.0), 2.0, -5.0, 0.0);
let e_high = Ellipse::from_radius_ratio(Point::new(0.0, 0.0), 2.0, 2.0, 0.0);
assert!(e_low.semi_major() > 0.0);
assert!(e_low.semi_minor() > 0.0);
assert!(e_high.semi_major() > 0.0);
assert!(e_high.semi_minor() > 0.0);
assert!((e_low.semi_minor() / e_low.semi_major() - 0.2).abs() < 0.01);
assert!((e_high.semi_minor() / e_high.semi_major() - 1.0).abs() < 1e-10);
}
#[test]
fn test_diagram_shape_to_params_is_geometric() {
use crate::geometry::traits::DiagramShape;
let e = Ellipse::new(Point::new(0.5, -1.5), 2.0, 0.5, 0.25);
let params = e.to_params();
assert_eq!(params.len(), 5);
assert_eq!(params[0], 0.5); assert_eq!(params[1], -1.5); assert_eq!(params[2], 2.0); assert_eq!(params[3], 0.5); assert_eq!(params[4], 0.25);
let back = Ellipse::from_params(¶ms);
assert_eq!(back.center(), e.center());
assert!((back.semi_major() - e.semi_major()).abs() < 1e-12);
assert!((back.semi_minor() - e.semi_minor()).abs() < 1e-12);
assert!((back.rotation() - e.rotation()).abs() < 1e-12);
}
#[test]
fn test_diagram_shape_optimizer_params_from_circle() {
use crate::geometry::traits::DiagramShape;
let params = Ellipse::optimizer_params_from_circle(1.0, 2.0, 3.0);
let log3 = 3.0_f64.ln();
assert_eq!(params.len(), 5);
assert_eq!(params[0], 1.0); assert_eq!(params[1], 2.0); assert!((params[2] - log3).abs() < 1e-12); assert!((params[3] - log3).abs() < 1e-12); assert_eq!(params[4], 0.0);
let e = Ellipse::from_optimizer_params(¶ms);
assert!((e.semi_major() - 3.0).abs() < 1e-10);
assert!((e.semi_minor() - 3.0).abs() < 1e-10);
}
#[test]
fn test_exclusive_regions_two_ellipses_partial_overlap() {
use crate::geometry::traits::DiagramShape;
let e1 = Ellipse::new(Point::new(0.0, 0.0), 4.0, 2.5, 0.0);
let e2 = Ellipse::new(Point::new(3.0, 1.0), 3.5, 2.0, PI / 6.0);
let areas = Ellipse::compute_exclusive_regions(&[e1, e2]);
let a_only = areas.get(&(1usize << 0)).copied().unwrap_or(0.0);
let b_only = areas.get(&(1usize << 1)).copied().unwrap_or(0.0);
let both = areas
.get(&((1usize << 0) | (1usize << 1)))
.copied()
.unwrap_or(0.0);
assert!(a_only >= 0.0 && b_only >= 0.0 && both >= 0.0);
let union = a_only + b_only + both;
assert!(union <= e1.area() + e2.area() + 1e-6);
let exact = e1.intersection_area(&e2);
let error = (both - exact).abs() / exact.max(1e-6);
assert!(
error < 0.1,
"Error too large: both={}, exact={}, error={:.2}%",
both,
exact,
error * 100.0
);
}
#[test]
fn test_exclusive_regions_two_ellipses_contained() {
use crate::geometry::traits::DiagramShape;
let outer = Ellipse::new(Point::new(0.0, 0.0), 6.0, 4.0, 0.2);
let inner = Ellipse::new(Point::new(0.5, 0.3), 2.0, 1.5, -0.3);
let areas = Ellipse::compute_exclusive_regions(&[outer, inner]);
let inner_only = areas.get(&(1usize << 1)).copied().unwrap_or(0.0);
let both = areas
.get(&((1usize << 0) | (1usize << 1)))
.copied()
.unwrap_or(0.0);
assert!((both - inner.area()).abs() < 1e-3);
assert!(inner_only < 1e-6);
}
#[test]
fn test_exclusive_regions_three_ellipses_basic() {
use crate::geometry::traits::DiagramShape;
let e1 = Ellipse::new(Point::new(-2.0, 0.0), 3.0, 2.0, 0.2);
let e2 = Ellipse::new(Point::new(2.0, 0.5), 3.2, 2.1, -0.1);
let e3 = Ellipse::new(Point::new(0.0, 1.5), 2.8, 1.8, 0.5);
let areas = Ellipse::compute_exclusive_regions(&[e1, e2, e3]);
let mask_all = (1usize << 0) | (1usize << 1) | (1usize << 2);
let all = areas.get(&mask_all).copied().unwrap_or(0.0);
let min_area = e1.area().min(e2.area()).min(e3.area());
assert!(all >= 0.0 && all <= min_area + 1e-6);
let union: f64 = areas.values().sum();
assert!(union <= e1.area() + e2.area() + e3.area() + 1e-3);
}
fn monte_carlo_pairwise_intersection(
e1: &Ellipse,
e2: &Ellipse,
n_samples: usize,
seed: u64,
) -> f64 {
use rand::rngs::StdRng;
use rand::{Rng, SeedableRng};
let mut rng = StdRng::seed_from_u64(seed);
let mut in_both = 0;
let bbox1 = e1.bounding_box();
let bbox2 = e2.bounding_box();
let (min1, max1) = bbox1.to_points();
let (min2, max2) = bbox2.to_points();
let x_min = min1.x().max(min2.x());
let x_max = max1.x().min(max2.x());
let y_min = min1.y().max(min2.y());
let y_max = max1.y().min(max2.y());
if x_min >= x_max || y_min >= y_max {
return 0.0; }
let bbox_area = (x_max - x_min) * (y_max - y_min);
for _ in 0..n_samples {
let x = x_min + (x_max - x_min) * rng.random::<f64>();
let y = y_min + (y_max - y_min) * rng.random::<f64>();
let p = Point::new(x, y);
if e1.contains_point(&p) && e2.contains_point(&p) {
in_both += 1;
}
}
(in_both as f64 / n_samples as f64) * bbox_area
}
#[test]
#[ignore = "slow stochastic validation"]
fn test_random_ellipse_intersections_monte_carlo() {
use rand::rngs::StdRng;
use rand::{Rng, SeedableRng};
const N_TESTS: usize = 100;
const MC_SAMPLES: usize = 50_000;
const MAX_RELATIVE_ERROR: f64 = 0.10;
let mut failures = Vec::new();
for seed in 0..N_TESTS {
let mut rng = StdRng::seed_from_u64(seed as u64);
let x1 = rng.random_range(-5.0..5.0);
let y1 = rng.random_range(-5.0..5.0);
let a1 = rng.random_range(0.5..3.0);
let b1 = rng.random_range(0.5..3.0);
let rot1 = rng.random_range(0.0..PI);
let x2 = x1 + rng.random_range(-2.0..2.0);
let y2 = y1 + rng.random_range(-2.0..2.0);
let a2 = rng.random_range(0.5..3.0);
let b2 = rng.random_range(0.5..3.0);
let rot2 = rng.random_range(0.0..PI);
let e1 = Ellipse::new(Point::new(x1, y1), a1, b1, rot1);
let e2 = Ellipse::new(Point::new(x2, y2), a2, b2, rot2);
let exact = e1.intersection_area(&e2);
let mc_estimate = monte_carlo_pairwise_intersection(&e1, &e2, MC_SAMPLES, seed as u64);
if exact < 0.01 && mc_estimate < 0.01 {
continue;
}
let max_val = exact.max(mc_estimate);
let error = (exact - mc_estimate).abs() / max_val;
if error > MAX_RELATIVE_ERROR {
let n_points = e1.intersection_points(&e2).len();
failures.push((seed, exact, mc_estimate, error, n_points));
}
}
assert!(
failures.is_empty(),
"Monte Carlo validation failed for {} of {} cases: {:?}",
failures.len(),
N_TESTS,
failures.iter().take(10).collect::<Vec<_>>()
);
}
#[test]
fn test_three_ellipse_complete_overlap() {
use crate::geometry::traits::DiagramShape;
let e1 = Ellipse::new(Point::new(0.0, 0.0), 0.564, 0.564, 0.0); let e2 = Ellipse::new(Point::new(0.0, 0.0), 0.564, 0.564, 0.0);
let e3 = Ellipse::new(Point::new(0.0, 0.0), 0.564, 0.564, 0.0);
let areas = Ellipse::compute_exclusive_regions(&[e1, e2, e3]);
let mask_all = (1usize << 0) | (1usize << 1) | (1usize << 2);
let all_three = areas.get(&mask_all).copied().unwrap_or(0.0);
assert!(
(all_three - 1.0).abs() < 0.01,
"Complete overlap should give area ~1.0, got {}",
all_three
);
}
#[test]
fn test_three_ellipse_partial_overlap_monte_carlo() {
use rand::rngs::StdRng;
use rand::{Rng, SeedableRng};
let seed = 42u64;
let mut rng = StdRng::seed_from_u64(seed);
let x1 = 0.0;
let y1 = 0.0;
let a1 = 1.5;
let b1 = 1.0;
let rot1 = 0.0;
let x2 = 1.0;
let y2 = 0.0;
let a2 = 1.5;
let b2 = 1.0;
let rot2 = PI / 6.0;
let x3 = 0.5;
let y3 = 1.0;
let a3 = 1.5;
let b3 = 1.0;
let rot3 = -PI / 6.0;
let e1 = Ellipse::new(Point::new(x1, y1), a1, b1, rot1);
let e2 = Ellipse::new(Point::new(x2, y2), a2, b2, rot2);
let e3 = Ellipse::new(Point::new(x3, y3), a3, b3, rot3);
let n_samples = 100_000;
let mut in_all_three = 0;
let bbox1 = e1.bounding_box();
let bbox2 = e2.bounding_box();
let bbox3 = e3.bounding_box();
let (min1, max1) = bbox1.to_points();
let (min2, max2) = bbox2.to_points();
let (min3, max3) = bbox3.to_points();
let x_min = min1.x().max(min2.x()).max(min3.x());
let x_max = max1.x().min(max2.x()).min(max3.x());
let y_min = min1.y().max(min2.y()).max(min3.y());
let y_max = max1.y().min(max2.y()).min(max3.y());
if x_min >= x_max || y_min >= y_max {
return;
}
let bbox_area = (x_max - x_min) * (y_max - y_min);
for _ in 0..n_samples {
let x = x_min + (x_max - x_min) * rng.random::<f64>();
let y = y_min + (y_max - y_min) * rng.random::<f64>();
let p = Point::new(x, y);
if e1.contains_point(&p) && e2.contains_point(&p) && e3.contains_point(&p) {
in_all_three += 1;
}
}
let mc_area = (in_all_three as f64 / n_samples as f64) * bbox_area;
use crate::geometry::traits::DiagramShape;
let areas = Ellipse::compute_exclusive_regions(&[e1, e2, e3]);
let mask_all = (1usize << 0) | (1usize << 1) | (1usize << 2);
let exact_area = areas.get(&mask_all).copied().unwrap_or(0.0);
if mc_area < 0.01 && exact_area < 0.01 {
return;
}
let error = (exact_area - mc_area).abs() / mc_area.max(exact_area);
assert!(
error < 0.15,
"Three-ellipse Monte Carlo error too large: {:.1}%",
error * 100.0
);
}
}