Crate grid1d

Crate grid1d 

Source
Expand description

One-dimensional grids and interval partitions with mathematical rigor and performance optimization.

This library provides a comprehensive framework for representing and manipulating one-dimensional grids and interval partitions with strong mathematical foundations, type safety, and performance optimizations. It serves as the backbone for numerical computations, finite element methods, adaptive mesh refinement, and scientific computing applications requiring robust 1D spatial discretizations.

§Core Design Philosophy

§Mathematical Correctness First

The library prioritizes mathematical rigor over convenience:

  • Type-safe intervals: Different interval types (IntervalClosed, IntervalOpen, etc.) encode boundary semantics at compile time
  • Domain validation: Runtime checks ensure grids completely partition their specified domains
  • Proper boundary handling: Open/closed boundary semantics are correctly preserved in all operations
  • Generic scalar support: Works with any scalar type implementing num_valid::RealScalar

§Multiple Floating-Point Backends via num-valid

grid1d leverages the num-valid library to provide multiple numerical backends with different precision and performance characteristics:

§Native Backend (Default)
§Arbitrary Precision Backend (Optional, activated with --features=rug)
  • Base types: rug::Float, rug::Complex from the GNU Multiple Precision Arithmetic Library
  • Validated types: num_valid::RealRugStrictFinite<N>
  • Performance: ⚡ Configurable, depends on precision setting
  • Precision: User-defined N-bit precision (e.g., 100, 256, 1024 bits)
  • Best for: Scientific computing requiring exact arithmetic, symbolic computation

All grid operations work seamlessly across backends - you can switch from standard to arbitrary precision by simply changing the scalar type, with no changes to your grid manipulation code.

use grid1d::{Grid1D, IntervalPartition, intervals::IntervalClosed, scalars::NumIntervals};
use num_valid::{RealNative64StrictFinite, RealScalar};
use try_create::TryNew;

// Standard precision grid
let domain = IntervalClosed::new(
    RealNative64StrictFinite::try_from_f64(0.0).unwrap(),
    RealNative64StrictFinite::try_from_f64(1.0).unwrap()
);
let standard_grid = Grid1D::uniform(domain, NumIntervals::try_new(100).unwrap());

use num_valid::RealRugStrictFinite;

// High precision grid (256-bit precision) - same API!
let hp_domain = IntervalClosed::new(
    RealRugStrictFinite::<256>::try_from_f64(0.0).unwrap(),
    RealRugStrictFinite::<256>::try_from_f64(1.0).unwrap()
);
let high_precision_grid = Grid1D::uniform(hp_domain, NumIntervals::try_new(100).unwrap());

// Both grids support identical operations
let point = RealRugStrictFinite::<256>::try_from_f64(0.5).unwrap();
let interval_id = high_precision_grid.find_interval_id_of_point(&point);

§Backend Selection Guide

BackendScalar TypePrecisionPerformanceUse Case
Nativef6453 bits⚡⚡⚡ MaximumProduction simulations, real-time
Native Validated (in Debug only)num_valid::RealNative64StrictFiniteInDebug53 bits⚡⚡⚡ Maximum (same as f64)Most applications
Native Validatednum_valid::RealNative64StrictFinite53 bits⚡⚡ Small penaltySafety-critical applications
Arbitrary Precisionnum_valid::RealRugStrictFinite<N>N bits⚡ ConfigurableScientific research, symbolic math

Enabling arbitrary precision: Add features = ["rug"] to your Cargo.toml:

[dependencies]
grid1d = { version = "0.1.0", features = ["rug"] } // change version to the latest available

§Performance Through Specialization

  • Uniform grids: O(1) point location with analytical formulas (Grid1DUniform)
  • Non-uniform grids: O(log n) binary search with optimized data structures (Grid1DNonUniform)
  • Memory efficiency: Zero-cost abstractions and cache-friendly memory layouts
  • SIMD-ready: Regular patterns enable vectorization optimizations

§Composability and Integration

  • Trait-based design: Core traits (IntervalPartition, HasDomain1D, HasCoords1D) enable generic programming
  • Unified interface: Grid1D enum provides seamless switching between uniform and non-uniform grids
  • Grid operations: Union, refinement, and intersection operations with bidirectional mappings

§Quick Start Guide

§Creating Grids

use grid1d::{
    Grid1D, HasCoords1D, IntervalPartition,
    intervals::{IntervalClosed, IntervalOpen},
    scalars::{NumIntervals, IntervalId},
};
use sorted_vec::partial::SortedSet;
use std::ops::Deref;
use try_create::TryNew;

// Create uniform grid with equal spacing
let domain = IntervalClosed::new(0.0, 1.0);
let uniform_grid = Grid1D::uniform(domain.clone(), NumIntervals::try_new(4).unwrap());
assert_eq!(uniform_grid.coords().deref(), &[0.0, 0.25, 0.5, 0.75, 1.0]);

// Create non-uniform grid with custom spacing
let coords = SortedSet::from_unsorted(vec![0.0, 0.1, 0.5, 0.9, 1.0]);
let non_uniform_grid = Grid1D::try_from_sorted(domain, coords).unwrap();

// Both grids implement the same IntervalPartition trait
assert_eq!(uniform_grid.num_intervals().as_ref(), &4);
assert_eq!(non_uniform_grid.num_intervals().as_ref(), &4);

§Point Location and Interval Access

use grid1d::{Grid1D, IntervalPartition, intervals::IntervalClosed, scalars::NumIntervals};
use try_create::TryNew;

let grid = Grid1D::uniform(IntervalClosed::new(0.0, 10.0), NumIntervals::try_new(10).unwrap());

// Efficient point location
let interval_id = grid.find_interval_id_of_point(&3.7);
println!("Point 3.7 is in interval {}", interval_id.as_ref());

// Access interval properties
let interval = grid.interval(&interval_id);
let length = grid.interval_length(&interval_id);
println!("Interval length: {}", length.as_ref());

// Batch point location
let points = vec![1.5, 4.2, 7.8, 9.1];
let intervals = grid.find_intervals_for_points(&points);

§Working with Different Scalar Types

use grid1d::{Grid1D, intervals::IntervalClosed, scalars::NumIntervals};
use num_valid::RealNative64StrictFiniteInDebug;
use try_create::TryNew;

// Use validated scalar types for debug safety
type Real = RealNative64StrictFiniteInDebug;
let domain = IntervalClosed::new(
    Real::try_new(0.0).unwrap(),
    Real::try_new(1.0).unwrap()
);
let grid = Grid1D::uniform(domain, NumIntervals::try_new(100).unwrap());

§Core Concepts

§Grid Types and Performance Characteristics

Grid TypePoint LocationMemoryBest For
Grid1DUniformO(1) analyticalO(n)Equal spacing, maximum performance
Grid1DNonUniformO(log n) binary searchO(n)Adaptive spacing, complex features
Grid1DUnionO(log n) on unified gridO(n+m)Multi-physics, grid combination

§Interval Partition Types

The library correctly handles different interval boundary semantics:

  • Closed intervals [a,b]: Include both endpoints
  • Open intervals (a,b): Exclude both endpoints
  • Semi-open intervals: [a,b) or (a,b] - include one endpoint
  • Automatic sub-interval construction: Preserves domain semantics in partitions

§Mathematical Properties and Guarantees

Every grid maintains these invariants:

  • Complete domain coverage: ⋃ intervals = domain (exact reconstruction)
  • Non-overlapping intervals: Intervals are disjoint except at boundaries
  • Unique point location: Every domain point belongs to exactly one interval
  • Boundary consistency: First/last coordinates match domain bounds exactly
  • Sorted coordinates: Points are in strict ascending order

§Advanced Features

§Adaptive Mesh Refinement

use grid1d::{*, intervals::*, scalars::*};
use try_create::TryNew;

let base_grid = Grid1D::uniform(IntervalClosed::new(0.0, 1.0), NumIntervals::try_new(4).unwrap());

// Uniform refinement: double resolution everywhere
let uniform_refinement = base_grid.refine_uniform(&PositiveNumPoints1D::try_new(1).unwrap());
assert_eq!(uniform_refinement.refined_grid().num_intervals().as_ref(), &8);

// Selective refinement: refine only specific intervals
let selective_plan = [
    (IntervalId::new(1), PositiveNumPoints1D::try_new(2).unwrap()), // 3 sub-intervals
    (IntervalId::new(3), PositiveNumPoints1D::try_new(1).unwrap()), // 2 sub-intervals
];
let selective_refinement = uniform_refinement.into_refined_grid().refine(&selective_plan);

§Grid Union for Multi-Physics

use grid1d::{*, intervals::*, scalars::*};
use sorted_vec::partial::SortedSet;
use try_create::TryNew;

let domain = IntervalClosed::new(0.0, 2.0);

// Physics A: coarse global grid
let grid_a = Grid1D::uniform(domain.clone(), NumIntervals::try_new(4).unwrap());

// Physics B: fine local grid
let coords_b = SortedSet::from_unsorted(vec![0.0, 0.3, 0.7, 1.1, 1.4, 2.0]);
let grid_b = Grid1D::try_from_sorted(domain, coords_b).unwrap();

// Create unified grid preserving mappings to both original grids
let union = Grid1DUnion::try_new(&grid_a, &grid_b).unwrap();
println!("Unified grid has {} intervals", union.num_refined_intervals().as_ref());

// Map data between grids
for (refined_id, a_id, b_id) in union.iter_interval_mappings() {
    println!("Unified interval {} maps to A[{}] and B[{}]",
             refined_id.as_ref(), a_id.as_ref(), b_id.as_ref());
}

§Domain-Specific Interval Semantics

use grid1d::{*, intervals::*, scalars::*};
use sorted_vec::partial::SortedSet;

// Closed domain: [0,2] includes both boundaries
let closed_grid = Grid1D::try_from_sorted(
    IntervalClosed::new(0.0, 2.0),
    SortedSet::from_unsorted(vec![0.0, 1.0, 2.0])
).unwrap();

// Open domain: (0,2) excludes both boundaries  
let open_grid = Grid1D::try_from_sorted(
    IntervalOpen::new(0.0, 2.0),
    SortedSet::from_unsorted(vec![0.0, 1.0, 2.0])
).unwrap();

// Different sub-interval types are automatically constructed
let closed_interval_0 = closed_grid.interval(&IntervalId::new(0)); // [0,1]
let open_interval_0 = open_grid.interval(&IntervalId::new(0));     // (0,1]

// Point containment reflects domain semantics
use grid1d::intervals::IntervalTrait;
// Note: This is a simplified example - actual containment checking requires pattern matching

§High-Performance Computations

use grid1d::{*, intervals::*, scalars::*};
use try_create::TryNew;

// Large uniform grid for finite difference methods
let grid = Grid1D::uniform(
    IntervalClosed::new(0.0, 1.0),
    NumIntervals::try_new(100_000).unwrap()
);

// O(1) point location for uniform grids
let start = std::time::Instant::now();
let test_points: Vec<f64> = (0..10_000).map(|i| i as f64 / 10_000.0).collect();
let intervals: Vec<IntervalId> = test_points.iter()
    .map(|&p| grid.find_interval_id_of_point(&p))
    .collect();
let elapsed = start.elapsed();
println!("Located {} points in {:?}", test_points.len(), elapsed);

§Performance Characteristics

§Time Complexity Summary

OperationUniform GridNon-Uniform GridGrid UnionNotes
Grid CreationO(n)O(n)O(n+m)Linear in point count
Point LocationO(1)O(log n)O(log n)Uniform grids use analytical formulas
Interval AccessO(1)O(1)O(1)Direct array indexing
IntersectionO(k)O(log n + k)O(log n + k)Where k = result size
Union Creation--O(n+m)Optimal merge algorithm
Grid RefinementO(n×r)O(n×r)O(n×r)Where r = refinement factor

§Space Complexity

  • Coordinates storage: O(n) for all grid types
  • Grid union mappings: O(n+m) additional space for bidirectional mappings
  • Refinement structures: O(refined intervals) for mapping preservation
  • Memory layout: Cache-friendly contiguous arrays for optimal performance

§Integration with other numerical Ecosystem

§Scalar Type Recommendations

Choose scalar types based on your performance and safety requirements:

use grid1d::{Grid1D, intervals::IntervalClosed, scalars::NumIntervals};
use num_valid::{RealNative64StrictFiniteInDebug, RealNative64StrictFinite};
use try_create::TryNew;

// Maximum performance (production)
let fast_grid = Grid1D::uniform(
    IntervalClosed::new(0.0_f64, 1.0_f64),
    NumIntervals::try_new(1000).unwrap()
);

// Balanced performance + debug safety (recommended)
type Real = RealNative64StrictFiniteInDebug;
let safe_grid = Grid1D::uniform(
    IntervalClosed::new(Real::try_new(0.0).unwrap(), Real::try_new(1.0).unwrap()),
    NumIntervals::try_new(1000).unwrap()
);

// Always validated (safety-critical)
type SafeReal = RealNative64StrictFinite;
let secure_grid = Grid1D::uniform(
    IntervalClosed::new(SafeReal::try_new(0.0).unwrap(), SafeReal::try_new(1.0).unwrap()),
    NumIntervals::try_new(1000).unwrap()
);

§Numerical Method Integration

use grid1d::{*, intervals::*, scalars::*};
use try_create::TryNew;

// Finite difference method setup
fn setup_finite_difference(domain: IntervalClosed<f64>, resolution: usize)
    -> (Grid1D<IntervalClosed<f64>>, Vec<f64>)
{
    let grid = Grid1D::uniform(domain, NumIntervals::try_new(resolution).unwrap());
    let solution = vec![0.0; grid.coords().len()];
    (grid, solution)
}

// Spectral method with periodic boundaries
fn setup_spectral_method(n_modes: usize) -> Grid1D<IntervalLowerClosedUpperOpen<f64>> {
    let domain = IntervalLowerClosedUpperOpen::new(0.0, 2.0 * std::f64::consts::PI);
    Grid1D::uniform(domain, NumIntervals::try_new(n_modes).unwrap())
}

§Error Handling and Validation

§Construction Error Types

The library provides detailed error information for debugging:

use grid1d::{*, intervals::*, scalars::*};
use sorted_vec::partial::SortedSet;

// Insufficient points
let coords = SortedSet::from_unsorted(vec![0.0]); // Only one point
match Grid1D::try_from_sorted(IntervalClosed::new(0.0, 1.0), coords) {
    Err(ErrorsGrid1D::RequiredAtLeastTwoDistinctPoints { num_points_provided, .. }) => {
        println!("Need at least 2 points, got {}", num_points_provided);
    }
    _ => unreachable!(),
}

// Boundary mismatches
let coords = SortedSet::from_unsorted(vec![0.1, 1.0]); // Wrong lower bound
match Grid1D::try_from_sorted(IntervalClosed::new(0.0, 1.0), coords) {
    Err(ErrorsGrid1D::LowerBoundMismatch { first_point, domain_lower_bound, .. }) => {
        println!("First point {} ≠ domain lower bound {}", first_point, domain_lower_bound);
    }
    _ => unreachable!(),
}

§Safe Construction Patterns

use grid1d::{*, intervals::*, scalars::*};
use sorted_vec::partial::SortedSet;

fn safe_grid_from_points(
    points: Vec<f64>,
    domain: IntervalClosed<f64>
) -> Result<Grid1D<IntervalClosed<f64>>, String> {
    if points.len() < 2 {
        return Err("Need at least 2 points".to_string());
    }
     
    let coords = SortedSet::from_unsorted(points);
    Grid1D::try_from_sorted(domain, coords)
        .map_err(|e| format!("Grid construction failed: {}", e))
}

§Best Practices and Recommendations

§When to Use Each Grid Type

Use Grid1DUniform when:

  • Finite difference methods requiring equal spacing
  • Spectral methods with uniform point distributions
  • Maximum performance is critical (O(1) point location)
  • Simple problems with uniform solution characteristics

Use Grid1DNonUniform when:

  • Adaptive mesh refinement based on solution features
  • Boundary layer resolution in fluid dynamics
  • Multi-scale problems requiring variable resolution
  • Complex geometries or irregular point distributions

Use Grid1DUnion when:

  • Combining grids from different physics solvers
  • Multi-physics coupling requiring consistent discretization
  • Grid overlay operations for interpolation
  • Preserving relationships between multiple grid levels

§Performance Optimization Tips

use grid1d::{*, intervals::*, scalars::*};
use try_create::TryNew;

// 1. Choose uniform grids when possible for O(1) point location
let uniform_grid = Grid1D::uniform(
    IntervalClosed::new(0.0, 1.0),
    NumIntervals::try_new(1000).unwrap()
);

// 2. Use batch operations for multiple points
let points = vec![0.1, 0.3, 0.7, 0.9];
let intervals = uniform_grid.find_intervals_for_points(&points); // More efficient than individual calls

// 3. Pre-allocate vectors for performance-critical loops
let mut results = Vec::with_capacity(points.len());
for &point in &points {
    results.push(uniform_grid.find_interval_id_of_point(&point));
}

// 4. Use direct coordinate access for sequential processing
let coords = uniform_grid.coords();
for i in 0..coords.len()-1 {
    let interval_length = coords[i+1] - coords[i];
    // Process interval...
}

§Memory Management

use grid1d::{*, intervals::*, scalars::*};
use try_create::TryNew;

// For large grids, consider the memory footprint
let large_grid = Grid1D::uniform(
    IntervalClosed::new(0.0, 1.0),
    NumIntervals::try_new(1_000_000).unwrap()
);

// Memory usage: ~8MB for coordinates (1M points × 8 bytes/f64)
println!("Grid memory: ~{} MB",
         large_grid.coords().len() * std::mem::size_of::<f64>() / 1_000_000);

// Use iterators to avoid creating intermediate vectors
let sum: f64 = large_grid.coords().iter().sum();
let mean = sum / large_grid.coords().len() as f64;

This module integrates with the other numerical ecosystem:

§Future Extensions

The design supports natural extensions to:

  • Multi-dimensional grids: 2D/3D generalizations using the same principles
  • Adaptive algorithms: Error-driven refinement strategies
  • Parallel processing: Thread-safe grid operations and domain decomposition
  • GPU acceleration: SIMD-friendly data layouts for vectorization
  • File I/O: Serialization for large-scale scientific computing

§Examples and Use Cases

§Computational Fluid Dynamics

use grid1d::{*, intervals::*, scalars::*};
use sorted_vec::partial::SortedSet;

// Create boundary layer mesh with fine spacing near walls
fn create_boundary_layer_mesh(domain_length: f64, boundary_thickness: f64)
    -> Grid1D<IntervalClosed<f64>>
{
    let mut points = vec![0.0];
     
    // Fine spacing in boundary layer
    for i in 1..=20 {
        let eta = (i as f64) / 20.0;
        let y = boundary_thickness * (eta * eta); // Quadratic clustering
        points.push(y);
    }
     
    // Coarse spacing in outer region
    for i in 1..=10 {
        let y = boundary_thickness + (domain_length - boundary_thickness) * (i as f64) / 10.0;
        points.push(y);
    }
     
    Grid1D::try_from_sorted(
        IntervalClosed::new(0.0, domain_length),
        SortedSet::from_unsorted(points)
    ).unwrap()
}

let cfd_grid = create_boundary_layer_mesh(1.0, 0.1);
println!("CFD grid: {} intervals, min spacing: {:.2e}",
         cfd_grid.num_intervals().as_ref(),
         cfd_grid.min_interval_length().as_ref());

§Adaptive Finite Element Analysis

use grid1d::{*, intervals::*, scalars::*};
use try_create::TryNew;

// Simulate adaptive refinement based on error indicators
fn adaptive_fem_cycle(
    mut grid: Grid1D<IntervalClosed<f64>>,
    error_threshold: f64,
    max_iterations: usize
) -> Grid1D<IntervalClosed<f64>> {
    for iteration in 0..max_iterations {
        // Compute error indicators (simplified)
        let errors: Vec<f64> = (0..*grid.num_intervals().as_ref())
            .map(|i| {
                // Simulate error based on interval size (larger intervals = higher error)
                grid.interval_length(&IntervalId::new(i)).as_ref() * 10.0
            })
            .collect();
         
        // Mark intervals for refinement
        let refinement_plan: Vec<(IntervalId, PositiveNumPoints1D)> = errors
            .iter()
            .enumerate()
            .filter_map(|(i, &error)| {
                if error > error_threshold {
                    Some((IntervalId::new(i), PositiveNumPoints1D::try_new(1).unwrap()))
                } else {
                    None
                }
            })
            .collect();
         
        if refinement_plan.is_empty() {
            println!("Converged after {} iterations", iteration);
            break;
        }
         
        // Apply refinement
        let refinement = grid.refine(&refinement_plan);
        grid = Grid1D::NonUniform(refinement.into_refined_grid());
         
        println!("Iteration {}: refined {} intervals, total: {}",
                 iteration, refinement_plan.len(), grid.num_intervals().as_ref());
    }
     
    grid
}

For comprehensive examples, API documentation, and advanced usage patterns, see the individual struct and trait documentation throughout this module.

Modules§

bounds
Mathematical bounds for 1D grid construction and interval arithmetic.
intervals
Comprehensive type-safe system for mathematical intervals with generic scalar support.
scalars
Type-safe scalar wrappers and numerical utilities for mathematical grid operations.

Structs§

Coords1D
Container for one-dimensional coordinates that are unique, sorted, and non-empty.
Coords1DInDomain
FINAL DOC A validated collection of one-dimensional coordinates guaranteed to lie within a specific interval domain.
Grid1DNonUniform
FINAL DOC A general-purpose, non-uniform partition of a finite interval with positive length.
Grid1DRefinement
A refinement of a 1D grid preserving bidirectional relationships between original and refined intervals.
Grid1DUniform
FINAL DOC A uniform partition of a finite interval with positive length into equal-sized sub-intervals.
Grid1DUnion
FINAL DOC Result of computing the union of two 1D grids with interval mappings.

Enums§

ErrorsCoords1D
Enum representing errors that can occur when constructing or manipulating one-dimensional coordinates.
ErrorsGeneralPoints1D
Enum representing errors that can occur when constructing or manipulating a collection of general 1D points.
ErrorsGrid1D
FINAL DOC Comprehensive error types for 1D grid construction and manipulation.
ErrorsGrid1DUnion
FINAL DOC Comprehensive error types for 1D grid union construction and operations.
Grid1D
An enum representing a one-dimensional grid, which defines a partition of a bounded interval with positive length.

Traits§

HasCoords1D
A trait for types that represent a set of one-dimensional coordinates.
HasDomain1D
Core trait for accessing one-dimensional domains.
IntervalInPartitionBuilder
A trait for building sub-intervals within a partition based on the domain type.
IntervalPartition
A trait for modeling a partition of an interval with finite, positive length as a set of adjacent non-overlapping intervals.
TransformedPoints1D
Trait representing a set of 1D points that are transformed between a parametric interval and a physical interval.

Functions§

linspace
Generates uniformly spaced points over a finite interval with positive length.
logspace
Generates logarithmically spaced points in the exponent domain.

Type Aliases§

Grid1DNonUniformRefinement
FINAL DOC A specialized refinement of a 1D grid where intervals are selectively subdivided with variable resolution.
Grid1DUniformRefinement
FINAL DOC A specialized refinement of a 1D grid where all intervals are uniformly subdivided.