grid1d 0.2.0

A mathematically rigorous, type-safe Rust library for 1D grid operations and interval partitions, supporting both native and arbitrary-precision numerics.
Documentation

grid1d

Crates.io Docs.rs License: MIT OR Apache-2.0 Pipeline Status Coverage GitLab last commit

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

grid1d is a comprehensive Rust library for representing and manipulating one-dimensional grids, interval partitions and interval arithmetic with strong mathematical foundations, type safety, and performance optimizations. Designed for numerical computations, finite element methods, adaptive mesh refinement, and scientific computing applications requiring robust 1D spatial discretizations.

๐ŸŽฏ Key Features

Mathematical Correctness First

  • Type-safe intervals: Different interval types ([a,b], (a,b), 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 correctly preserved in all operations
  • Generic scalar support: Works with any scalar type implementing num_valid::RealScalar

Multiple Floating-Point Backends

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

Backend Scalar Type Precision Performance Best For
Native f64 53 bits โšกโšกโšก Maximum Production simulations, real-time
Native Validated (Debug) RealNative64StrictFiniteInDebug 53 bits โšกโšกโšก Same as f64 Recommended for most applications
Native Validated RealNative64StrictFinite 53 bits โšกโšก Small overhead Safety-critical applications
Arbitrary Precision RealRugStrictFinite<N> N bits โšก Configurable Scientific research, symbolic math

Performance Through Specialization

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

Advanced Grid Operations

  • Adaptive mesh refinement: Uniform and selective refinement with bidirectional mappings
  • Grid union: Combine grids from different physics for multi-physics simulations
  • Interval arithmetic: Complete set of interval operations with proper boundary handling
  • Type-safe scalars: Prevent common numerical bugs with validated scalar types

๐Ÿš€ Quick Start

Prerequisites

This library requires Rust nightly for advanced language features:

# rust-toolchain.toml (automatically configured)
[toolchain]
channel = "nightly"

Why nightly is required:

  • Advanced error handling (#![feature(error_generic_member_access)]): Better debugging and error context
  • Future optimizations: Enables cutting-edge Rust developments for numerical computing

Installation

Add this to your Cargo.toml:

[dependencies]
grid1d = "0.1.2"  # Change to the latest version
try_create = "0.1.2"  # Required for TryNew trait
sorted-vec = "0.8"  # Required for non-uniform grids

# For arbitrary precision arithmetic (optional)
grid1d = { version = "0.1.2", features = ["rug"] }
num-valid = { version = "0.2.1", features = ["rug"] }

Your First Grid: Complete Tutorial

Let's build a complete example from scratch, step by step:

Step 1: Create a Simple Uniform Grid

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

fn main() {
    // Define the domain: [0, 10]
    let domain = IntervalClosed::new(0.0, 10.0);
    
    // Create a grid with 10 equal intervals
    let grid = Grid1D::uniform(domain, NumIntervals::try_new(10).unwrap());
    
    // Grid has 11 points: 0, 1, 2, ..., 10
    println!("Grid has {} points", grid.coords().len());
}

Step 2: Locate Points in the Grid

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

fn main() {
    let domain = IntervalClosed::new(0.0, 10.0);
    let grid = Grid1D::uniform(domain, NumIntervals::try_new(10).unwrap());
    
    // Find which interval contains a point
    let point = 3.7;
    let interval_id = grid.find_interval_id_of_point(&point);
    println!("Point {} is in interval {}", point, interval_id.as_ref());
    
    // Get the actual interval
    let interval = grid.interval(&interval_id);
    println!("Interval: {:?}", interval);
    
    // Get interval length
    let length = grid.interval_length(&interval_id);
    println!("Interval length: {}", length.as_ref());
}

Step 3: Work with Non-Uniform Grids

use grid1d::{Grid1D, intervals::IntervalClosed};
use sorted_vec::partial::SortedSet;

fn main() {
    // Custom point distribution (e.g., clustered near x=0)
    let coords = SortedSet::from_unsorted(vec![
        0.0, 0.1, 0.2, 0.5, 1.0, 2.0, 5.0, 10.0
    ]);
    
    let domain = IntervalClosed::new(0.0, 10.0);
    let grid = Grid1D::try_from_sorted(domain, coords).unwrap();
    
    println!("Non-uniform grid with {} intervals", grid.num_intervals().as_ref());
}

Step 4: Choose the Right Scalar Type

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

fn main() {
    // Use validated types for safety (recommended)
    type Real = RealNative64StrictFiniteInDebug;
    
    let domain = IntervalClosed::new(
        Real::try_new(0.0).unwrap(),
        Real::try_new(10.0).unwrap()
    );
    
    let grid = Grid1D::uniform(domain, NumIntervals::try_new(10).unwrap());
    
    // This will catch NaN/Infinity in debug builds
    let safe_point = Real::try_new(3.7).unwrap();
    let interval_id = grid.find_interval_id_of_point(&safe_point);
}

Decision Guide: Which Grid Type Should I Use?

// โœ… Use Grid1DUniform when:
// - All intervals have equal spacing
// - Maximum performance is critical (O(1) point location)
// - Regular discretization (e.g., finite differences on uniform mesh)

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

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

// โœ… Use Grid1DNonUniform when:
// - Adaptive spacing needed (e.g., boundary layers, shock capturing)
// - Point clustering in regions of interest
// - Arbitrary point distributions

use sorted_vec::partial::SortedSet;

let clustered_coords = SortedSet::from_unsorted(vec![/* custom distribution */]);
let non_uniform = Grid1D::try_from_sorted(
    IntervalClosed::new(0.0, 1.0),
    clustered_coords
).unwrap();

// โœ… Use Grid1DUnion when:
// - Combining grids from multiple physics or solvers
// - Multi-physics simulations with different discretizations
// - Need to maintain mappings between original and unified grids

let grid_a = Grid1D::uniform(/* ... */);
let grid_b = Grid1D::uniform(/* ... */);
// let union = Grid1DUnion::try_new(&grid_a, &grid_b).unwrap();

Decision Guide: Which Scalar Type Should I Use?

Your Need Recommended Type Example
Maximum performance, trusted inputs f64 Production code with validated inputs
โญ Recommended default RealNative64StrictFiniteInDebug Development + production (same performance as f64)
Safety-critical applications RealNative64StrictFinite Medical, aerospace, financial
Arbitrary precision needed RealRugStrictFinite<N> Scientific computing, symbolic math
// Maximum performance
let fast_domain = IntervalClosed::new(0.0_f64, 1.0_f64);

// Recommended: Safe in debug, fast in release
use num_valid::RealNative64StrictFiniteInDebug;
type Real = RealNative64StrictFiniteInDebug;
let safe_domain = IntervalClosed::new(
    Real::try_new(0.0).unwrap(),
    Real::try_new(1.0).unwrap()
);

// Always validated
use num_valid::RealNative64StrictFinite;
type SafeReal = RealNative64StrictFinite;
let validated_domain = IntervalClosed::new(
    SafeReal::try_new(0.0).unwrap(),
    SafeReal::try_new(1.0).unwrap()
);

Basic Usage

use grid1d::{
    Grid1D, IntervalPartition, HasCoords1D,
    intervals::IntervalClosed,
    scalars::{NumIntervals, IntervalId},
};
use sorted_vec::partial::SortedSet;
use std::ops::Deref;  // Needed for .deref()
use try_create::TryNew;  // Needed for .try_new()

// 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);

// Efficient point location
let interval_id = uniform_grid.find_interval_id_of_point(&0.3);
let interval = uniform_grid.interval(&interval_id);

Working with Different Scalar Types

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

// Performance-optimal with debug safety (recommended)
type OptimalReal = RealNative64StrictFiniteInDebug;
let optimal_grid = Grid1D::uniform(
    IntervalClosed::new(
        OptimalReal::try_new(0.0).unwrap(),
        OptimalReal::try_new(1.0).unwrap()
    ),
    NumIntervals::try_new(4).unwrap()
);

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

Arbitrary Precision with Rug

For applications requiring arbitrary precision arithmetic:

// Add to Cargo.toml:
// grid1d = { version = "0.1.0", features = ["rug"] }

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

// 256-bit precision arithmetic
type HighPrecisionReal = RealRugStrictFinite<256>;

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

// All grid operations work seamlessly with arbitrary precision
let hp_1 = HighPrecisionReal::try_new(1.0).unwrap();
let hp_3 = HighPrecisionReal::try_new(3.0).unwrap();
let point = hp_1 / hp_3; // 1. / 3.
let interval_id = high_precision_grid.find_interval_id_of_point(&point);

๐Ÿ”ฌ Advanced Features

Adaptive Mesh Refinement

use grid1d::{*, intervals::*, scalars::*};
use std::collections::BTreeMap;  // Needed for refinement plan
use try_create::TryNew;  // Needed for .try_new()

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 = BTreeMap::from([
    (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;  // Needed for .try_new()

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());
}

Working with Different Interval Types

use grid1d::intervals::*;

let closed = IntervalClosed::new(0.0, 1.0);         // [0, 1]
let open = IntervalOpen::new(0.0, 1.0);             // (0, 1)
let half_open = IntervalLowerClosedUpperOpen::new(0.0, 1.0); // [0, 1)
let unbounded = IntervalLowerClosedUpperUnbounded::new(0.0); // [0, +โˆž)

// All implement IntervalTrait for generic operations
assert!(closed.contains_point(&0.0));
assert!(!open.contains_point(&0.0));

๐Ÿ—๏ธ Core Concepts

Grid Types and Performance Characteristics

Grid Type Point Location Memory Best For
Grid1DUniform O(1) analytical O(n) Equal spacing, maximum performance
Grid1DNonUniform O(log n) binary search O(n) Adaptive spacing, complex features
Grid1DUnion O(log n) on unified grid O(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

๐Ÿ“– Examples

Finite Difference Methods

use grid1d::{*, intervals::*, scalars::*};
use std::ops::Deref;  // Needed for accessing coords
use try_create::TryNew;  // Needed for .try_new()

let grid = Grid1D::uniform(IntervalClosed::new(0.0, 1.0), NumIntervals::try_new(100).unwrap());
let coords = grid.coords().deref();
let dx = coords[1] - coords[0];

// Initialize solution vector
let mut solution = vec![0.0; coords.len()];
// ... (set initial/boundary conditions)

// Apply finite difference stencils (interior points only)
for i in 1..coords.len()-1 {
    let d2u_dx2 = (solution[i-1] - 2.0*solution[i] + solution[i+1]) / (dx * dx);
    // Time stepping or iteration logic...
    // e.g., solution_new[i] = solution[i] + dt * d2u_dx2;
}

Spectral Methods

use grid1d::{*, intervals::*, scalars::*};
use try_create::TryNew;  // Needed for .try_new()
use std::ops::Deref;  // Needed for accessing coords

// Periodic domain for Fourier spectral methods
let domain = IntervalLowerClosedUpperOpen::new(0.0, 2.0 * std::f64::consts::PI);
let grid = Grid1D::uniform(domain, NumIntervals::try_new(256).unwrap());

// Perfect for FFT algorithms
let dx = grid.coords()[1] - grid.coords()[0];
assert!((dx - 2.0 * std::f64::consts::PI / 256.0).abs() < 1e-15);

Boundary Layer Refinement

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", cfd_grid.num_intervals().as_ref());

๐Ÿ› ๏ธ Development

Building

# Standard build (uses nightly automatically)
cargo build

# With arbitrary precision support
cargo build --features=rug

# Run tests
cargo test

# Run tests with arbitrary precision
cargo test --features=rug

# Generate documentation
cargo doc --open

Requirements

  • Rust: Nightly toolchain (specified in rust-toolchain.toml)

Core Dependencies

  • [num-valid]: Validated scalar arithmetic
  • [sorted-vec]: Efficient sorted container implementation
  • [serde]: Serialization support for grid structures
  • [nutype]: Type-safe wrapper generation
  • [derive-more]: Ergonomic derive macros
  • [getset]: Automatic getter generation

Optional Dependencies

  • [rug]: Arbitrary precision arithmetic via the rug library
    • Enable with --features=rug
    • Provides num_valid::RealRugStrictFinite<N> types for N-bit precision

๐Ÿงช Testing

The library includes comprehensive test suites:

# Run all tests
cargo test

# Run tests with arbitrary precision
cargo test --features=rug

# Run property-based tests
cargo test --test property_tests

# Run benchmarks
cargo bench

Test Coverage

  • Unit tests: Core functionality for all grid types
  • Integration tests: Cross-module interactions
  • Property-based tests: Mathematical invariants and edge cases
  • Benchmarks: Performance regression detection
  • Documentation tests: All code examples in docs are tested

Development Setup

  1. Clone the repository:

    git clone git@gitlab.com:max.martinelli/grid1d.git
    cd grid1d
    
  2. Ensure you have the required dependencies:

    • Rust nightly toolchain
  3. Run the test suite:

    cargo test --all-features
    
  4. Check formatting and linting:

    cargo fmt --check
    cargo clippy --all-features
    

Contribution Guidelines

  • Mathematical correctness: All changes must preserve mathematical guarantees
  • Performance: Performance-critical paths must maintain or improve benchmarks
  • Documentation: All public APIs must be thoroughly documented with examples
  • Testing: New features require comprehensive test coverage
  • Type safety: Leverage Rust's type system for compile-time correctness

๐Ÿ“„ License

This project is licensed under the MIT License - see the LICENSE file for details.

๐Ÿ”— Related Projects

  • [num-valid]: Generic validated scalar arithmetic
  • [ndarray]: N-dimensional arrays for Rust
  • [nalgebra]: Linear algebra library
  • [faer]: High-performance linear algebra
  • [polars]: Fast DataFrame library

๐Ÿ“š Documentation

For comprehensive API documentation, mathematical foundations, and advanced usage patterns:

Quick Links

๐ŸŽฏ Use Cases

grid1d is designed for:

Scientific Computing

  • Computational Fluid Dynamics: Boundary layer meshes, shock capturing
  • Finite Element Analysis: Adaptive mesh refinement, domain decomposition
  • Spectral Methods: Fourier transforms, Chebyshev polynomials
  • Numerical Integration: Quadrature rules, adaptive integration

Engineering Applications

  • Heat Transfer: Non-uniform grids for boundary layers
  • Structural Analysis: Stress concentration regions
  • Signal Processing: Non-uniform sampling, filter design
  • Control Systems: State-space discretization

Data Science

  • Time Series: Irregular time grids, event-driven sampling
  • Interpolation: Scattered data, adaptive sampling
  • Numerical Methods: Root finding, optimization
  • Machine Learning: Feature engineering, data preprocessing

License

Copyright 2023-2025, C.N.R. - Consiglio Nazionale delle Ricerche

Licensed under either of

at your option.

Unless you explicitly state otherwise, any contribution intentionally submitted for inclusion in this project by you, as defined in the Apache-2.0 license, shall be dual licensed as above, without any additional terms or conditions.

License Notice for Optional Feature Dependencies (LGPL-3.0 Compliance)

If you enable the rug feature, this project will depend on the rug library, which is licensed under the LGPL-3.0 license. Activating this feature may introduce LGPL-3.0 requirements to your project. Please review the terms of the LGPL-3.0 license to ensure compliance.