IntervalFinitePositiveLength

Enum IntervalFinitePositiveLength 

Source
pub enum IntervalFinitePositiveLength<RealType: RealScalar> {
    Closed(IntervalClosed<RealType>),
    Open(IntervalOpen<RealType>),
    LowerOpenUpperClosed(IntervalLowerOpenUpperClosed<RealType>),
    LowerClosedUpperOpen(IntervalLowerClosedUpperOpen<RealType>),
}
Expand description

FINAL DOC Container for all intervals with finite, positive, measurable length.

The IntervalFinitePositiveLength enum represents intervals where the upper bound is strictly greater than the lower bound, resulting in a positive measure and uncountably infinite cardinality. This enum encompasses all bounded intervals with non-zero length, distinguished by their boundary inclusion behavior.

§Mathematical Foundation

All variants represent intervals of the form {x ∈ ℝ : a < x < b} or similar, where the specific boundary conditions depend on the variant. The fundamental mathematical properties shared by all variants are:

  • Positive Measure: μ(I) = |b - a| > 0 for all variants
  • Non-empty Interior: All intervals contain infinitely many interior points
  • Boundedness: Both endpoints are finite real numbers
  • Connectedness: All intervals are connected sets with no gaps
  • Uncountable Cardinality: Each interval contains uncountably many points

§Generic Over Any num_valid::RealScalar Type

The IntervalFinitePositiveLength enum works seamlessly with any scalar type implementing num_valid::RealScalar:

Scalar TypePerformanceValidationBest For
f64⚡⚡⚡ Maximum❌ NoneTrusted input, maximum speed
RealNative64StrictFiniteInDebug⚡⚡⚡ Same as f64Debug onlyRecommended for most uses
RealNative64StrictFinite⚡⚡ Small overhead✅ AlwaysSafety-critical applications
RealRugStrictFinite⚡ Precision-dependent✅ AlwaysArbitrary precision needs (available from the num-valid crate when compiled with --features=rug)

§Enum Variants and Boundary Semantics

The four variants differ only in their boundary inclusion behavior:

§IntervalFinitePositiveLength::Closed - [a, b]

Both endpoints included

  • Mathematical notation: [a, b] = {x ∈ ℝ : a ≤ x ≤ b}
  • Boundary behavior: Contains both a and b
  • Topology: Compact set (closed and bounded)
  • Supremum/Infimum: Both achieved at the endpoints
  • Use cases: Probability ranges, integration bounds, physical constraints

§IntervalFinitePositiveLength::Open - (a, b)

Both endpoints excluded

  • Mathematical notation: (a, b) = {x ∈ ℝ : a < x < b}
  • Boundary behavior: Excludes both a and b
  • Topology: Open set, not compact
  • Supremum/Infimum: Exist but not achieved
  • Use cases: Strict inequalities, limit domains, convergence regions

§IntervalFinitePositiveLength::LowerOpenUpperClosed - (a, b]

Lower endpoint excluded, upper endpoint included

  • Mathematical notation: (a, b] = {x ∈ ℝ : a < x ≤ b}
  • Boundary behavior: Excludes a, includes b
  • Topology: Neither open nor closed (half-open)
  • Use cases: Right-continuous functions, cumulative distributions

§IntervalFinitePositiveLength::LowerClosedUpperOpen - [a, b)

Lower endpoint included, upper endpoint excluded

  • Mathematical notation: [a, b) = {x ∈ ℝ : a ≤ x < b}
  • Boundary behavior: Includes a, excludes b
  • Topology: Neither open nor closed (half-open)
  • Use cases: Left-continuous functions, array indexing, partitions

§Construction Patterns

use grid1d::intervals::*;
use num_valid::RealNative64StrictFiniteInDebug;
use try_create::TryNew;

// Use the recommended performance-optimal type
type Real = RealNative64StrictFiniteInDebug;

// Create specific interval types, then convert
let closed = IntervalClosed::new(
    Real::try_new(0.0).unwrap(),
    Real::try_new(1.0).unwrap()
);
let positive: IntervalFinitePositiveLength<Real> = closed.into();

let open = IntervalOpen::new(
    Real::try_new(-1.0).unwrap(),
    Real::try_new(1.0).unwrap()
);
let positive: IntervalFinitePositiveLength<Real> = open.into();

§Direct Construction (Advanced)

use grid1d::intervals::*;
use try_create::New;

// Build variants directly
let closed = IntervalFinitePositiveLength::Closed(
    IntervalClosed::new(0.0, 1.0)
);

let open = IntervalFinitePositiveLength::Open(
    IntervalOpen::new(0.0, 1.0)
);

let half_open_left = IntervalFinitePositiveLength::LowerOpenUpperClosed(
    IntervalLowerOpenUpperClosed::new(0.0, 1.0)
);

let half_open_right = IntervalFinitePositiveLength::LowerClosedUpperOpen(
    IntervalLowerClosedUpperOpen::new(0.0, 1.0)
);

§Core Operations

All IntervalFinitePositiveLength instances support the complete IntervalTrait interface plus specialized operations for positive-length intervals:

§Point Containment Testing

use grid1d::intervals::*;
use try_create::New;

let closed = IntervalFinitePositiveLength::Closed(
    IntervalClosed::new(0.0, 1.0)
);
let open = IntervalFinitePositiveLength::Open(
    IntervalOpen::new(0.0, 1.0)
);

// Closed interval includes endpoints
assert!(closed.contains_point(&0.0));    // Lower bound included
assert!(closed.contains_point(&0.5));    // Interior point
assert!(closed.contains_point(&1.0));    // Upper bound included
assert!(!closed.contains_point(&2.0));   // Outside interval

// Open interval excludes endpoints
assert!(!open.contains_point(&0.0));     // Lower bound excluded
assert!(open.contains_point(&0.5));      // Interior point
assert!(!open.contains_point(&1.0));     // Upper bound excluded
assert!(!open.contains_point(&2.0));     // Outside interval

§Geometric Properties

use grid1d::intervals::*;
use try_create::New;

let interval = IntervalFinitePositiveLength::Closed(
    IntervalClosed::new(-2.0, 3.0)
);

// All variants have identical geometric properties
assert_eq!(interval.length().into_inner(), 5.0);           // |b - a| = |3 - (-2)| = 5
assert_eq!(interval.mid_point(), 0.5);        // (a + b) / 2 = (-2 + 3) / 2 = 0.5
assert!(!interval.is_symmetric());            // a + b ≠ 0

// Bounds access
assert_eq!(interval.lower_bound_value(), &-2.0);
assert_eq!(interval.upper_bound_value(), &3.0);

// Consume to get bounds
let (lower, upper) = interval.into_bounds_pair();
assert_eq!((lower, upper), (-2.0, 3.0));

§Interval Containment Testing

use grid1d::intervals::*;
use try_create::New;

let outer_closed = IntervalFinitePositiveLength::Closed(
    IntervalClosed::new(0.0, 2.0)
);
let inner_open = IntervalFinitePositiveLength::Open(
    IntervalOpen::new(0.5, 1.5)
);
let boundary_half_open = IntervalFinitePositiveLength::LowerClosedUpperOpen(
    IntervalLowerClosedUpperOpen::new(0.0, 1.0)
);

// Closed intervals can contain open intervals
assert!(outer_closed.contains_interval(&inner_open));           // [0,2] ⊇ (0.5,1.5) ✓
assert!(outer_closed.contains_interval(&boundary_half_open));   // [0,2] ⊇ [0,1) ✓

// Open intervals cannot contain intervals that include excluded boundaries
assert!(!inner_open.contains_interval(&outer_closed));          // (0.5,1.5) ⊉ [0,2] ✗
assert!(!inner_open.contains_interval(&boundary_half_open));    // (0.5,1.5) ⊉ [0,1) ✗

§Intersection Operations

use grid1d::intervals::*;

let closed = IntervalFinitePositiveLength::Closed(
    IntervalClosed::new(0.0, 2.0)
);
let open = IntervalFinitePositiveLength::Open(
    IntervalOpen::new(1.0, 3.0)
);

if let Some(intersection) = closed.intersection(&open) {
    // Result: (1.0, 2.0] - most restrictive boundary combination
    match intersection {
        Interval::FiniteLength(IntervalFiniteLength::PositiveLength(positive)) => {
            println!("Intersection length: {}", positive.length()); // 1.0
             
            // Can extract specific type if needed
            if let IntervalFinitePositiveLength::LowerOpenUpperClosed(half_open) = positive {
                assert_eq!(half_open.lower_bound_value(), &1.0);
                assert_eq!(half_open.upper_bound_value(), &2.0);
            }
        }
        _ => unreachable!("Positive-length intervals always have positive-length intersections"),
    }
}

// Edge case: intersection at single point becomes singleton
let touching_left = IntervalFinitePositiveLength::Closed(
    IntervalClosed::new(0.0, 1.0)
);
let touching_right = IntervalFinitePositiveLength::Closed(
    IntervalClosed::new(1.0, 2.0)
);

if let Some(intersection) = touching_left.intersection(&touching_right) {
    // Result: singleton at point 1.0
    match intersection {
        Interval::FiniteLength(IntervalFiniteLength::ZeroLength(singleton)) => {
            assert_eq!(singleton.value(), &1.0);
        }
        _ => panic!("Expected singleton intersection"),
    }
}

§Pattern Matching and Variant Analysis

§Distinguishing Between Boundary Types

use grid1d::intervals::*;
use try_create::New;

fn analyze_boundary_behavior(interval: &IntervalFinitePositiveLength<f64>) -> String {
    let lower = *interval.lower_bound_value();
    let upper = *interval.upper_bound_value();
     
    match interval {
        IntervalFinitePositiveLength::Closed(_) => {
            format!("[{}, {}] - both endpoints included", lower, upper)
        }
        IntervalFinitePositiveLength::Open(_) => {
            format!("({}, {}) - both endpoints excluded", lower, upper)
        }
        IntervalFinitePositiveLength::LowerOpenUpperClosed(_) => {
            format!("({}, {}] - lower open, upper closed", lower, upper)
        }
        IntervalFinitePositiveLength::LowerClosedUpperOpen(_) => {
            format!("[{}, {}) - lower closed, upper open", lower, upper)
        }
    }
}

let closed = IntervalFinitePositiveLength::Closed(IntervalClosed::new(1.0, 2.0));
let open = IntervalFinitePositiveLength::Open(IntervalOpen::new(1.0, 2.0));
let half_open_left = IntervalFinitePositiveLength::LowerOpenUpperClosed(
    IntervalLowerOpenUpperClosed::new(1.0, 2.0)
);
let half_open_right = IntervalFinitePositiveLength::LowerClosedUpperOpen(
    IntervalLowerClosedUpperOpen::new(1.0, 2.0)
);

println!("{}", analyze_boundary_behavior(&closed));       // "[1, 2] - both endpoints included"
println!("{}", analyze_boundary_behavior(&open));         // "(1, 2) - both endpoints excluded"
println!("{}", analyze_boundary_behavior(&half_open_left)); // "(1, 2] - lower open, upper closed"
println!("{}", analyze_boundary_behavior(&half_open_right)); // "[1, 2) - lower closed, upper open"

§Boundary Endpoint Testing

use grid1d::intervals::*;

fn test_boundary_inclusion(interval: &IntervalFinitePositiveLength<f64>) -> (bool, bool) {
    let lower = *interval.lower_bound_value();
    let upper = *interval.upper_bound_value();
     
    let includes_lower = interval.contains_point(&lower);
    let includes_upper = interval.contains_point(&upper);
     
    (includes_lower, includes_upper)
}

let closed = IntervalFinitePositiveLength::Closed(IntervalClosed::new(0.0, 1.0));
let open = IntervalFinitePositiveLength::Open(IntervalOpen::new(0.0, 1.0));
let left_open = IntervalFinitePositiveLength::LowerOpenUpperClosed(
    IntervalLowerOpenUpperClosed::new(0.0, 1.0)
);
let right_open = IntervalFinitePositiveLength::LowerClosedUpperOpen(
    IntervalLowerClosedUpperOpen::new(0.0, 1.0)
);

assert_eq!(test_boundary_inclusion(&closed), (true, true));     // [0,1]: both included
assert_eq!(test_boundary_inclusion(&open), (false, false));     // (0,1): both excluded
assert_eq!(test_boundary_inclusion(&left_open), (false, true)); // (0,1]: lower excluded, upper included
assert_eq!(test_boundary_inclusion(&right_open), (true, false)); // [0,1): lower included, upper excluded

§Mathematical Applications

§Integration Domains

use grid1d::intervals::*;
use num_valid::RealScalar;

/// Compute definite integral using midpoint rule (simplified example)
fn integrate_midpoint<T: RealScalar + Clone>(
    f: impl Fn(T) -> T,
    domain: &IntervalFinitePositiveLength<T>,
    num_subdivisions: usize
) -> T {
    let length = domain.length().into_inner();
    let step = length.clone() / T::try_from_f64(num_subdivisions as f64).unwrap();
    let mut sum = T::zero();
     
    for i in 0..num_subdivisions {
        let offset = T::try_from_f64(i as f64).unwrap() * step.clone() + step.clone() / T::try_from_f64(2.0).unwrap();
        let x = domain.lower_bound_value().clone() + offset;
        sum = sum + f(x);
    }
     
    sum * step
}

// Works with any boundary type - the integral is the same
let closed_domain = IntervalFinitePositiveLength::Closed(IntervalClosed::new(0.0, 1.0));
let open_domain = IntervalFinitePositiveLength::Open(IntervalOpen::new(0.0, 1.0));

let f = |x: f64| x * x; // f(x) = x²

let integral_closed = integrate_midpoint(&f, &closed_domain, 1000);
let integral_open = integrate_midpoint(&f, &open_domain, 1000);

// Both should approximate ∫₀¹ x² dx = 1/3 ≈ 0.333...
assert!((integral_closed - 1.0/3.0).abs() < 0.01);
assert!((integral_open - 1.0/3.0).abs() < 0.01);

§Probability Theory and Measure

use grid1d::intervals::*;

/// Represents a continuous uniform distribution
struct UniformDistribution {
    support: IntervalFinitePositiveLength<f64>
}

impl UniformDistribution {
    fn new(support: IntervalFinitePositiveLength<f64>) -> Self {
        Self { support }
    }
     
    /// Probability density function
    fn pdf(&self, x: f64) -> f64 {
        if self.support.contains_point(&x) {
            1.0 / self.support.length().into_inner() // Uniform density
        } else {
            0.0 // Outside support
        }
    }
     
    /// Cumulative distribution function  
    fn cdf(&self, x: f64) -> f64 {
        let lower = *self.support.lower_bound_value();
        let upper = *self.support.upper_bound_value();
         
        if x <= lower {
            0.0
        } else if x >= upper {
            1.0
        } else {
            (x - lower) / self.support.length().into_inner()
        }
    }
     
    /// Sample from the distribution (simplified)
    fn sample(&self) -> f64 {
        let lower = *self.support.lower_bound_value();
        let upper = *self.support.upper_bound_value();
         
        // In real implementation, use proper random number generator
        lower + (upper - lower) * 0.5 // Deterministic "sample" for testing
    }
}

// Works with any interval type
let uniform_closed = UniformDistribution::new(
    IntervalFinitePositiveLength::Closed(IntervalClosed::new(0.0, 2.0))
);
let uniform_open = UniformDistribution::new(
    IntervalFinitePositiveLength::Open(IntervalOpen::new(0.0, 2.0))
);

// Same density and behavior regardless of boundary type
assert_eq!(uniform_closed.pdf(1.0), 0.5);  // 1 / (2 - 0) = 0.5
assert_eq!(uniform_open.pdf(1.0), 0.5);    // Same density

assert_eq!(uniform_closed.cdf(1.0), 0.5);  // (1 - 0) / (2 - 0) = 0.5
assert_eq!(uniform_open.cdf(1.0), 0.5);    // Same CDF for interior points

§Optimization and Root Finding

use grid1d::intervals::*;

/// Simple bisection method for root finding
fn bisection_method(
    f: impl Fn(f64) -> f64,
    mut interval: IntervalFinitePositiveLength<f64>,
    tolerance: f64
) -> Option<f64> {
    let mut iterations = 0;
    const MAX_ITERATIONS: usize = 100;
     
    // Check that function changes sign over interval
    let f_lower = f(*interval.lower_bound_value());
    let f_upper = f(*interval.upper_bound_value());
     
    if f_lower * f_upper > 0.0 {
        return None; // No sign change, no root guaranteed
    }
     
    while interval.length().into_inner() > tolerance && iterations < MAX_ITERATIONS {
        let midpoint = interval.mid_point();
        let f_mid = f(midpoint);
         
        if f_mid.abs() < tolerance {
            return Some(midpoint);
        }
         
        let lower = *interval.lower_bound_value();
        let upper = *interval.upper_bound_value();
         
        // Update interval to half containing the root
        if f_lower * f_mid < 0.0 {
            // Root in left half
            interval = IntervalFinitePositiveLength::Closed(
                IntervalClosed::new(lower, midpoint)
            );
        } else {
            // Root in right half  
            interval = IntervalFinitePositiveLength::Closed(
                IntervalClosed::new(midpoint, upper)
            );
        }
         
        iterations += 1;
    }
     
    Some(interval.mid_point())
}

// Find root of f(x) = x² - 2 (should be √2 ≈ 1.414)
let search_interval = IntervalFinitePositiveLength::Closed(
    IntervalClosed::new(1.0, 2.0)
);

let root = bisection_method(|x| x * x - 2.0, search_interval, 1e-6).unwrap();
assert!((root - 2.0_f64.sqrt()).abs() < 1e-6);

§Type Conversion and Specialization

§Converting to Specific Interval Types

use grid1d::intervals::*;

let positive = IntervalFinitePositiveLength::Closed(
    IntervalClosed::new(0.0, 1.0)
);

// Convert to general interval enums
let finite: IntervalFiniteLength<f64> = positive.clone().into();
let general: Interval<f64> = positive.clone().into();

// Extract specific types when needed
match positive {
    IntervalFinitePositiveLength::Closed(closed_interval) => {
        println!("Working with closed interval: {:?}", closed_interval);
        // Access closed-interval-specific behavior if needed
        let (lower, upper) = closed_interval.into_bounds_pair();
        assert_eq!((lower, upper), (0.0, 1.0));
    }
    IntervalFinitePositiveLength::Open(open_interval) => {
        println!("Working with open interval: {:?}", open_interval);
    }
    IntervalFinitePositiveLength::LowerOpenUpperClosed(half_open) => {
        println!("Working with left half-open interval: {:?}", half_open);
    }
    IntervalFinitePositiveLength::LowerClosedUpperOpen(half_open) => {
        println!("Working with right half-open interval: {:?}", half_open);
    }
}

§Extracting Bounds

use grid1d::intervals::*;

fn extract_bounds_any_variant(
    interval: IntervalFinitePositiveLength<f64>
) -> (f64, f64) {
    // Works uniformly across all variants
    interval.into_bounds_pair()
}

let closed = IntervalFinitePositiveLength::Closed(IntervalClosed::new(1.0, 2.0));
let open = IntervalFinitePositiveLength::Open(IntervalOpen::new(3.0, 4.0));
let half_open = IntervalFinitePositiveLength::LowerOpenUpperClosed(
    IntervalLowerOpenUpperClosed::new(5.0, 6.0)
);

assert_eq!(extract_bounds_any_variant(closed), (1.0, 2.0));
assert_eq!(extract_bounds_any_variant(open), (3.0, 4.0));
assert_eq!(extract_bounds_any_variant(half_open), (5.0, 6.0));

§Performance Characteristics

§Memory Layout

  • Storage: Two scalar values plus boundary type discriminant
  • Size: Typically 16-24 bytes for f64, depending on alignment
  • Alignment: Natural alignment of the scalar type
  • Enum overhead: Small discriminant tag (1-8 bytes)

§Operation Complexity

OperationTime ComplexitySpace ComplexityNotes
Point queriesO(1)O(1)2 comparisons maximum
Containment testsO(1)O(1)Boundary logic only
IntersectionsO(1)O(1)Pure computation, may allocate result
Geometric calculationsO(1)O(1)Length, midpoint, etc.
Pattern matchingO(1)O(1)Compile-time optimized
Bound extractionO(1)O(1)Zero-cost access

§Scalar Type Performance Impact

use grid1d::intervals::*;
use num_valid::{RealNative64StrictFiniteInDebug, RealNative64StrictFinite, RealScalar};

// These have identical performance in release builds:
type FastPositiveInterval = IntervalFinitePositiveLength<f64>;
type OptimalPositiveInterval = IntervalFinitePositiveLength<RealNative64StrictFiniteInDebug>;

// This has small validation overhead:
type SafePositiveInterval = IntervalFinitePositiveLength<RealNative64StrictFinite>;

// Performance test function (generic over scalar type)
fn benchmark_positive_operations<T: RealScalar + Clone>(
    interval1: &IntervalFinitePositiveLength<T>,
    interval2: &IntervalFinitePositiveLength<T>
) {
    // All operations compile to identical assembly for f64 and RealNative64StrictFiniteInDebug
    let _contains = interval1.contains_interval(interval2);
    let _intersection = interval1.intersection(interval2);
    let _length = interval1.length().into_inner();        // O(1) subtraction
    let _midpoint = interval1.mid_point();   // O(1) arithmetic
    let _symmetric = interval1.is_symmetric(); // O(1) comparison
     
    // Pattern matching is zero-cost
    match interval1 {
        IntervalFinitePositiveLength::Closed(_) => { /* ... */ }
        IntervalFinitePositiveLength::Open(_) => { /* ... */ }
        IntervalFinitePositiveLength::LowerOpenUpperClosed(_) => { /* ... */ }
        IntervalFinitePositiveLength::LowerClosedUpperOpen(_) => { /* ... */ }
    }
}

§Integration with Library Ecosystem

§With Numerical Integration

use grid1d::intervals::*;

/// Adaptive quadrature integration
fn adaptive_integrate(
    f: impl Fn(f64) -> f64 + Copy,
    domain: &IntervalFinitePositiveLength<f64>,
    tolerance: f64
) -> f64 {
    // The integration algorithm is identical regardless of boundary type
    // because the measure is the same and boundary points have measure zero
     
    let length = domain.length().into_inner();
    let midpoint = domain.mid_point();
     
    // Simpson's rule approximation (simplified)
    let lower = *domain.lower_bound_value();
    let upper = *domain.upper_bound_value();
     
    (length / 6.0) * (f(lower) + 4.0 * f(midpoint) + f(upper))
}

// All boundary types give the same integral (up to numerical precision)
let closed_domain = IntervalFinitePositiveLength::Closed(
    IntervalClosed::new(0.0, 1.0)
);
let open_domain = IntervalFinitePositiveLength::Open(
    IntervalOpen::new(0.0, 1.0)
);

let f = |x: f64| x * x;
let integral_closed = adaptive_integrate(f, &closed_domain, 1e-6);
let integral_open = adaptive_integrate(f, &open_domain, 1e-6);

// Should be approximately equal (boundary points don't affect integral)
assert!((integral_closed - integral_open).abs() < 1e-10);

§With Mesh Generation

use grid1d::intervals::*;

/// Generate uniform mesh points within any positive-length interval
fn generate_uniform_mesh(
    domain: &IntervalFinitePositiveLength<f64>,
    num_points: usize
) -> Vec<f64> {
    if num_points <= 1 {
        return vec![domain.mid_point()];
    }
     
    let lower = *domain.lower_bound_value();
    let upper = *domain.upper_bound_value();
    let step = (upper - lower) / (num_points - 1) as f64;
     
    (0..num_points)
        .map(|i| lower + i as f64 * step)
        .collect()
}

/// Generate mesh respecting boundary conditions
fn generate_boundary_aware_mesh(
    domain: &IntervalFinitePositiveLength<f64>,
    num_interior_points: usize
) -> Vec<f64> {
    let lower = *domain.lower_bound_value();
    let upper = *domain.upper_bound_value();
     
    let mut points = Vec::new();
     
    // Include boundary points based on interval type
    match domain {
        IntervalFinitePositiveLength::Closed(_) => {
            points.push(lower);
        }
        IntervalFinitePositiveLength::LowerClosedUpperOpen(_) => {
            points.push(lower);
        }
        _ => {} // Open boundaries not included
    }
     
    // Add interior points
    for i in 1..=num_interior_points {
        let t = i as f64 / (num_interior_points + 1) as f64;
        points.push(lower + t * (upper - lower));
    }
     
    // Include upper boundary if appropriate
    match domain {
        IntervalFinitePositiveLength::Closed(_) => {
            points.push(upper);
        }
        IntervalFinitePositiveLength::LowerOpenUpperClosed(_) => {
            points.push(upper);
        }
        _ => {} // Open boundaries not included
    }
     
    points
}

§Best Practices and Recommendations

§When to Use IntervalFinitePositiveLength

  • Numerical analysis: When you need positive-length domains but want to abstract boundary behavior
  • Generic algorithms: Writing code that works with any finite positive interval
  • API design: Accepting any bounded interval while excluding singletons
  • Mathematical modeling: When the boundary type may vary but the domain must be non-degenerate

§When to Use Specific Variants

  • Closed: Integration bounds, probability ranges, physical constraints with inclusive limits
  • Open: Limit domains, convergence regions, strict inequality constraints
  • Half-open: Partitions, cumulative functions, discrete-continuous interfaces
use grid1d::intervals::*;
use into_inner::IntoInner;
use num_valid::{RealNative64StrictFiniteInDebug, functions::Rounding};
use num::Zero;

// Define type aliases for your domain
type Real = RealNative64StrictFiniteInDebug;
type PositiveInterval = IntervalFinitePositiveLength<Real>;
type ClosedInterval = IntervalClosed<Real>;

// Pattern: Handle all variants uniformly when possible
fn compute_interval_statistics(interval: &PositiveInterval) -> (Real, Real, Real) {
    let length = interval.length().into_inner();
    let midpoint = interval.mid_point();
    let lower = interval.lower_bound_value().clone();
     
    (length, midpoint, lower) // Same calculation for all variants
}

// Pattern: Specialize behavior only when boundary semantics matter
fn count_integer_endpoints(interval: &PositiveInterval) -> usize {
    let lower = interval.lower_bound_value().clone();
    let upper = interval.upper_bound_value().clone();
    let lower_int = lower.clone().kernel_floor().into_inner() as i64;
    let upper_int = upper.clone().kernel_ceil().into_inner() as i64;
     
    let mut count = 0;
     
    // Check lower boundary
    if lower == lower_int as f64 && interval.contains_point(&lower) {
        count += 1;
    }
     
    // Check upper boundary  
    if upper == upper_int as f64 && interval.contains_point(&upper) {
        count += 1;
    }
     
    count
}

// Pattern: Convert to most specific type when needed
fn require_closed_interval(interval: PositiveInterval) -> Result<ClosedInterval, String> {
    match interval {
        IntervalFinitePositiveLength::Closed(closed) => Ok(closed),
        _ => Err("Algorithm requires closed interval for boundary conditions".to_string()),
    }
}

§Error Handling and Edge Cases

§Construction Validation

use grid1d::intervals::*;

// All variants validate that lower < upper during construction
assert!(IntervalClosed::try_new(1.0, 0.0).is_err());  // Invalid bounds
assert!(IntervalOpen::try_new(1.0, 1.0).is_err());    // Zero length
assert!(IntervalLowerOpenUpperClosed::try_new(2.0, 1.0).is_err()); // Inverted

// In debug mode, invalid construction panics
#[cfg(debug_assertions)]
{
    std::panic::catch_unwind(|| {
        IntervalClosed::new(1.0, 0.0); // Panics in debug
    }).expect_err("Should panic on invalid bounds");
}

§Numerical Precision Considerations

use grid1d::intervals::*;

// CAUTION: Very small intervals near floating-point precision limits
let tiny_interval = IntervalFinitePositiveLength::Closed(
    IntervalClosed::new(1.0, 1.0 + f64::EPSILON)
);

// Length is positive but very small
assert!(tiny_interval.length().into_inner() > 0.0);
assert!(tiny_interval.length().into_inner() < 1e-15);

// Midpoint calculation may have precision issues
let midpoint = tiny_interval.mid_point();
let expected_mid = 1.0 + f64::EPSILON / 2.0;

// Use appropriate tolerance for comparisons
assert!((midpoint - expected_mid).abs() < f64::EPSILON);

§Common Pitfalls and Solutions

use grid1d::intervals::*;

// WRONG: Assuming all intervals include their endpoints
fn naive_endpoint_processing(interval: &IntervalFinitePositiveLength<f64>) -> (f64, f64) {
    let lower = *interval.lower_bound_value();
    let upper = *interval.upper_bound_value();
     
    // This ignores boundary semantics!
    (lower, upper) // May include points not actually in the interval
}

// CORRECT: Check boundary inclusion when it matters
fn proper_endpoint_processing(interval: &IntervalFinitePositiveLength<f64>) -> Vec<f64> {
    let lower = *interval.lower_bound_value();
    let upper = *interval.upper_bound_value();
     
    let mut endpoints = Vec::new();
     
    if interval.contains_point(&lower) {
        endpoints.push(lower);
    }
    if interval.contains_point(&upper) {
        endpoints.push(upper);
    }
     
    endpoints
}

let closed = IntervalFinitePositiveLength::Closed(IntervalClosed::new(0.0, 1.0));
let open = IntervalFinitePositiveLength::Open(IntervalOpen::new(0.0, 1.0));

assert_eq!(proper_endpoint_processing(&closed), vec![0.0, 1.0]); // Both endpoints
assert!(proper_endpoint_processing(&open).is_empty()); // No endpoints

§Mathematical Correctness Guarantees

The IntervalFinitePositiveLength enum maintains all essential mathematical properties:

  • Positive measure: μ(I) = |upper - lower| > 0 for all variants
  • Containment transitivity: If A ⊆ B and B ⊆ C, then A ⊆ C
  • Intersection commutativity: A ∩ B = B ∩ A
  • Intersection associativity: (A ∩ B) ∩ C = A ∩ (B ∩ C)
  • Monotonicity: If A ⊆ B, then μ(A) ≤ μ(B)
  • Boundary consistency: Containment and intersection operations respect boundary semantics
  • Geometric invariance: Length, midpoint, and symmetry properties are boundary-independent

These properties are preserved regardless of the scalar type used and are guaranteed by the underlying implementations of the IntervalTrait and IntervalFinitePositiveLengthTrait.

Variants§

§

Closed(IntervalClosed<RealType>)

A closed interval [a,b] where both endpoints are included.

This variant represents the most restrictive form of inclusion, where both boundary points a and b are considered part of the interval. The interval contains all real numbers x such that a ≤ x ≤ b.

Mathematical properties:

  • Compactness: Always compact (closed and bounded)
  • Boundary inclusion: Both a and b are contained
  • Supremum/Infimum: Both achieved at the endpoints
  • Typical use cases: Integration bounds, probability ranges, physical constraints
§

Open(IntervalOpen<RealType>)

An open interval (a,b) where both endpoints are excluded.

This variant represents the most permissive form of exclusion, where neither boundary point a nor b is considered part of the interval. The interval contains all real numbers x such that a < x < b.

Mathematical properties:

  • Openness: Open set in the real line topology
  • Boundary exclusion: Neither a nor b are contained
  • Supremum/Infimum: Exist but are not achieved
  • Typical use cases: Limit domains, convergence regions, strict inequalities
§

LowerOpenUpperClosed(IntervalLowerOpenUpperClosed<RealType>)

A left half-open interval (a,b] where the lower endpoint is excluded and the upper endpoint is included.

This variant represents a mixed boundary condition where the lower bound a is not part of the interval but the upper bound b is included. The interval contains all real numbers x such that a < x ≤ b.

Mathematical properties:

  • Boundary behavior: Excludes a, includes b
  • Topology: Neither open nor closed (half-open)
  • Right-continuity: Natural for right-continuous functions
  • Typical use cases: Cumulative distributions, partition elements
§

LowerClosedUpperOpen(IntervalLowerClosedUpperOpen<RealType>)

A right half-open interval [a,b) where the lower endpoint is included and the upper endpoint is excluded.

This variant represents a mixed boundary condition where the lower bound a is part of the interval but the upper bound b is not included. The interval contains all real numbers x such that a ≤ x < b.

Mathematical properties:

  • Boundary behavior: Includes a, excludes b
  • Topology: Neither open nor closed (half-open)
  • Left-continuity: Natural for left-continuous functions
  • Typical use cases: Array indexing ranges, interval partitions, discrete-continuous boundaries

Trait Implementations§

Source§

impl<RealType: Clone + RealScalar> Clone for IntervalFinitePositiveLength<RealType>

Source§

fn clone(&self) -> IntervalFinitePositiveLength<RealType>

Returns a duplicate of the value. Read more
1.0.0 · Source§

fn clone_from(&mut self, source: &Self)

Performs copy-assignment from source. Read more
Source§

impl<RealType: RealScalar> Contains<RealType> for IntervalFinitePositiveLength<RealType>

Source§

fn contains_point(&self, x: &RealType) -> bool

Returns true if the 1D point/coordinate x is contained in the current IntervalFinitePositiveLength.

Source§

fn contains_interval<IntervalType: IntervalTrait<RealType = RealType>>( &self, other: &IntervalType, ) -> bool

Tests whether this interval completely contains another interval. Read more
Source§

impl<RealType: Debug + RealScalar> Debug for IntervalFinitePositiveLength<RealType>

Source§

fn fmt(&self, f: &mut Formatter<'_>) -> Result

Formats the value using the given formatter. Read more
Source§

impl<'de, RealType> Deserialize<'de> for IntervalFinitePositiveLength<RealType>
where RealType: for<'a> Deserialize<'a> + RealScalar,

Source§

fn deserialize<__D>(__deserializer: __D) -> Result<Self, __D::Error>
where __D: Deserializer<'de>,

Deserialize this value from the given Serde deserializer. Read more
Source§

impl<RealType: RealScalar> From<IntervalBounded<RealType, Closed, Closed>> for IntervalFinitePositiveLength<RealType>

Source§

fn from(interval: IntervalClosed<RealType>) -> Self

Converts to this type from the input type.
Source§

impl<RealType: RealScalar> From<IntervalBounded<RealType, Closed, Open>> for IntervalFinitePositiveLength<RealType>

Source§

fn from(interval: IntervalLowerClosedUpperOpen<RealType>) -> Self

Converts to this type from the input type.
Source§

impl<RealType: RealScalar> From<IntervalBounded<RealType, Open, Closed>> for IntervalFinitePositiveLength<RealType>

Source§

fn from(interval: IntervalLowerOpenUpperClosed<RealType>) -> Self

Converts to this type from the input type.
Source§

impl<RealType: RealScalar> From<IntervalBounded<RealType, Open, Open>> for IntervalFinitePositiveLength<RealType>

Source§

fn from(interval: IntervalOpen<RealType>) -> Self

Converts to this type from the input type.
Source§

impl<RealType: RealScalar> From<IntervalFinitePositiveLength<RealType>> for Interval<RealType>

Source§

fn from(interval: IntervalFinitePositiveLength<RealType>) -> Self

Converts to this type from the input type.
Source§

impl<RealType: RealScalar> From<IntervalFinitePositiveLength<RealType>> for IntervalFiniteLength<RealType>

Source§

fn from(interval: IntervalFinitePositiveLength<RealType>) -> Self

Converts to this type from the input type.
Source§

impl<RealType: RealScalar> From<SubIntervalInPartition<RealType>> for IntervalFinitePositiveLength<RealType>

Source§

fn from(src: SubIntervalInPartition<RealType>) -> Self

Converts to this type from the input type.
Source§

impl<RealType: RealScalar> HasLowerBound for IntervalFinitePositiveLength<RealType>

Source§

type LowerBoundValue = RealType

Source§

fn lower_bound_value(&self) -> &RealType

Return the lower bound of the interval. Read more
Source§

impl<RealType: RealScalar> HasUpperBound for IntervalFinitePositiveLength<RealType>

Source§

type UpperBoundValue = RealType

Source§

fn upper_bound_value(&self) -> &RealType

Return the upper bound of the interval. Read more
Source§

impl<RealType: RealScalar> IntervalFinitePositiveLengthTrait for IntervalFinitePositiveLength<RealType>

Source§

fn into_bounds_pair(self) -> (Self::RealType, Self::RealType)

Return the inner bounds of the interval as a pair in which the first entry of the pair is the lower bound and the second entry of the pair is the upper bound.

§Note

This function consumes Self (no cloning of the bounds is needed).

Source§

fn length(&self) -> PositiveRealScalar<Self::RealType>

Return the length of the interval. Read more
Source§

fn mid_point(&self) -> Self::RealType

Return the midpoint of the interval. Read more
Source§

fn is_symmetric(&self) -> bool

Return true if the interval is symmetric about the origin. Read more
Source§

fn interior(self) -> IntervalOpen<Self::RealType>

Convert this interval to an open interval with the same bounds. Read more
Source§

fn closure(self) -> IntervalClosed<Self::RealType>

Convert this interval to a closed interval with the same bounds. Read more
Source§

impl<RealType: RealScalar> IntervalInPartitionBuilder for IntervalFinitePositiveLength<RealType>

Source§

fn build_first_interval( &self, lower_bound: RealType, upper_bound: RealType, ) -> IntervalFinitePositiveLength<RealType>

Constructs the first sub-interval in a partition. Read more
Source§

fn build_last_interval( &self, lower_bound: RealType, upper_bound: RealType, ) -> IntervalFinitePositiveLength<RealType>

Constructs the last sub-interval in a partition. Read more
Source§

impl<RealType: RealScalar> IntervalOperations for IntervalFinitePositiveLength<RealType>

Source§

type RealType = RealType

The scalar type used for interval bounds. Read more
Source§

fn intersection<IntervalType: IntervalTrait<RealType = Self::RealType>>( &self, other: &IntervalType, ) -> Option<Interval<Self::RealType>>

Computes the intersection of this interval with another interval. Read more
Source§

fn union<IntervalType: IntervalTrait<RealType = Self::RealType>>( &self, other: &IntervalType, ) -> IntervalUnion<Self::RealType>
where Interval<Self::RealType>: From<Self> + From<IntervalType>,

Computes the union of this interval with another. Read more
Source§

impl<RealType: PartialEq + RealScalar> PartialEq for IntervalFinitePositiveLength<RealType>

Source§

fn eq(&self, other: &IntervalFinitePositiveLength<RealType>) -> bool

Tests for self and other values to be equal, and is used by ==.
1.0.0 · Source§

fn ne(&self, other: &Rhs) -> bool

Tests for !=. The default implementation is almost always sufficient, and should not be overridden without very good reason.
Source§

impl<RealType> Serialize for IntervalFinitePositiveLength<RealType>
where RealType: Serialize + RealScalar,

Source§

fn serialize<__S>(&self, __serializer: __S) -> Result<__S::Ok, __S::Error>
where __S: Serializer,

Serialize this value into the given Serde serializer. Read more
Source§

impl<RealType: RealScalar> TryFrom<Interval<RealType>> for IntervalFinitePositiveLength<RealType>

Source§

type Error = ErrorsIntervalConversion<RealType>

The type returned in the event of a conversion error.
Source§

fn try_from(interval: Interval<RealType>) -> Result<Self, Self::Error>

Performs the conversion.
Source§

impl<RealType: RealScalar> TryFrom<IntervalFinitePositiveLength<RealType>> for IntervalClosed<RealType>

Source§

type Error = ErrorsIntervalConversion<RealType>

The type returned in the event of a conversion error.
Source§

fn try_from( interval_finite_positive_length: IntervalFinitePositiveLength<RealType>, ) -> Result<Self, Self::Error>

Performs the conversion.
Source§

impl<RealType: RealScalar> TryFrom<IntervalFinitePositiveLength<RealType>> for IntervalLowerClosedUpperOpen<RealType>

Source§

type Error = ErrorsIntervalConversion<RealType>

The type returned in the event of a conversion error.
Source§

fn try_from( interval_finite_positive_length: IntervalFinitePositiveLength<RealType>, ) -> Result<Self, Self::Error>

Performs the conversion.
Source§

impl<RealType: RealScalar> TryFrom<IntervalFinitePositiveLength<RealType>> for IntervalLowerOpenUpperClosed<RealType>

Source§

type Error = ErrorsIntervalConversion<RealType>

The type returned in the event of a conversion error.
Source§

fn try_from( interval_finite_positive_length: IntervalFinitePositiveLength<RealType>, ) -> Result<Self, Self::Error>

Performs the conversion.
Source§

impl<RealType: RealScalar> TryFrom<IntervalFinitePositiveLength<RealType>> for IntervalOpen<RealType>

Source§

type Error = ErrorsIntervalConversion<RealType>

The type returned in the event of a conversion error.
Source§

fn try_from( interval_finite_positive_length: IntervalFinitePositiveLength<RealType>, ) -> Result<Self, Self::Error>

Performs the conversion.
Source§

impl<RealType: Eq + RealScalar> Eq for IntervalFinitePositiveLength<RealType>

Source§

impl<RealType: RealScalar> IntervalTrait for IntervalFinitePositiveLength<RealType>

Source§

impl<RealType: RealScalar> StructuralPartialEq for IntervalFinitePositiveLength<RealType>

Auto Trait Implementations§

§

impl<RealType> Freeze for IntervalFinitePositiveLength<RealType>
where <RealType as RealScalar>::RawReal: Sized, RealType: Freeze,

§

impl<RealType> RefUnwindSafe for IntervalFinitePositiveLength<RealType>
where <RealType as RealScalar>::RawReal: Sized, RealType: RefUnwindSafe,

§

impl<RealType> Send for IntervalFinitePositiveLength<RealType>
where <RealType as RealScalar>::RawReal: Sized,

§

impl<RealType> Sync for IntervalFinitePositiveLength<RealType>
where <RealType as RealScalar>::RawReal: Sized,

§

impl<RealType> Unpin for IntervalFinitePositiveLength<RealType>
where <RealType as RealScalar>::RawReal: Sized, RealType: Unpin,

§

impl<RealType> UnwindSafe for IntervalFinitePositiveLength<RealType>
where <RealType as RealScalar>::RawReal: Sized, RealType: UnwindSafe,

Blanket Implementations§

Source§

impl<T> Any for T
where T: 'static + ?Sized,

Source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
Source§

impl<T> Az for T

Source§

fn az<Dst>(self) -> Dst
where T: Cast<Dst>,

Casts the value.
Source§

impl<T> Borrow<T> for T
where T: ?Sized,

Source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
Source§

impl<T> BorrowMut<T> for T
where T: ?Sized,

Source§

fn borrow_mut(&mut self) -> &mut T

Mutably borrows from an owned value. Read more
Source§

impl<Src, Dst> CastFrom<Src> for Dst
where Src: Cast<Dst>,

Source§

fn cast_from(src: Src) -> Dst

Casts the value.
Source§

impl<T> CheckedAs for T

Source§

fn checked_as<Dst>(self) -> Option<Dst>
where T: CheckedCast<Dst>,

Casts the value.
Source§

impl<Src, Dst> CheckedCastFrom<Src> for Dst
where Src: CheckedCast<Dst>,

Source§

fn checked_cast_from(src: Src) -> Option<Dst>

Casts the value.
Source§

impl<T> CloneToUninit for T
where T: Clone,

Source§

unsafe fn clone_to_uninit(&self, dest: *mut u8)

🔬This is a nightly-only experimental API. (clone_to_uninit)
Performs copy-assignment from self to dest. Read more
Source§

impl<T> From<T> for T

Source§

fn from(t: T) -> T

Returns the argument unchanged.

Source§

impl<T, U> Into<U> for T
where U: From<T>,

Source§

fn into(self) -> U

Calls U::from(self).

That is, this conversion is whatever the implementation of From<T> for U chooses to do.

Source§

impl<T> OverflowingAs for T

Source§

fn overflowing_as<Dst>(self) -> (Dst, bool)
where T: OverflowingCast<Dst>,

Casts the value.
Source§

impl<Src, Dst> OverflowingCastFrom<Src> for Dst
where Src: OverflowingCast<Dst>,

Source§

fn overflowing_cast_from(src: Src) -> (Dst, bool)

Casts the value.
Source§

impl<T> SaturatingAs for T

Source§

fn saturating_as<Dst>(self) -> Dst
where T: SaturatingCast<Dst>,

Casts the value.
Source§

impl<Src, Dst> SaturatingCastFrom<Src> for Dst
where Src: SaturatingCast<Dst>,

Source§

fn saturating_cast_from(src: Src) -> Dst

Casts the value.
Source§

impl<T> ToOwned for T
where T: Clone,

Source§

type Owned = T

The resulting type after obtaining ownership.
Source§

fn to_owned(&self) -> T

Creates owned data from borrowed data, usually by cloning. Read more
Source§

fn clone_into(&self, target: &mut T)

Uses borrowed data to replace owned data, usually by cloning. Read more
Source§

impl<T, U> TryFrom<U> for T
where U: Into<T>,

Source§

type Error = Infallible

The type returned in the event of a conversion error.
Source§

fn try_from(value: U) -> Result<T, <T as TryFrom<U>>::Error>

Performs the conversion.
Source§

impl<T, U> TryInto<U> for T
where U: TryFrom<T>,

Source§

type Error = <U as TryFrom<T>>::Error

The type returned in the event of a conversion error.
Source§

fn try_into(self) -> Result<U, <U as TryFrom<T>>::Error>

Performs the conversion.
Source§

impl<T> UnwrappedAs for T

Source§

fn unwrapped_as<Dst>(self) -> Dst
where T: UnwrappedCast<Dst>,

Casts the value.
Source§

impl<Src, Dst> UnwrappedCastFrom<Src> for Dst
where Src: UnwrappedCast<Dst>,

Source§

fn unwrapped_cast_from(src: Src) -> Dst

Casts the value.
Source§

impl<V, T> VZip<V> for T
where V: MultiLane<T>,

Source§

fn vzip(self) -> V

Source§

impl<T> WrappingAs for T

Source§

fn wrapping_as<Dst>(self) -> Dst
where T: WrappingCast<Dst>,

Casts the value.
Source§

impl<Src, Dst> WrappingCastFrom<Src> for Dst
where Src: WrappingCast<Dst>,

Source§

fn wrapping_cast_from(src: Src) -> Dst

Casts the value.
Source§

impl<T> DeserializeOwned for T
where T: for<'de> Deserialize<'de>,