grid1d
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)
[]
= "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:
[]
= "0.1.2" # Change to the latest version
= "0.1.2" # Required for TryNew trait
= "0.8" # Required for non-uniform grids
# For arbitrary precision arithmetic (optional)
= { = "0.1.2", = ["rug"] }
= { = "0.2.1", = ["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 ;
use TryNew;
Step 2: Locate Points in the Grid
use ;
use TryNew;
Step 3: Work with Non-Uniform Grids
use ;
use SortedSet;
Step 4: Choose the Right Scalar Type
use ;
use RealNative64StrictFiniteInDebug;
use TryNew;
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 ;
use TryNew;
let uniform = uniform;
// โ
Use Grid1DNonUniform when:
// - Adaptive spacing needed (e.g., boundary layers, shock capturing)
// - Point clustering in regions of interest
// - Arbitrary point distributions
use SortedSet;
let clustered_coords = from_unsorted;
let non_uniform = try_from_sorted.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 = uniform;
let grid_b = 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 = new;
// Recommended: Safe in debug, fast in release
use RealNative64StrictFiniteInDebug;
type Real = RealNative64StrictFiniteInDebug;
let safe_domain = new;
// Always validated
use RealNative64StrictFinite;
type SafeReal = RealNative64StrictFinite;
let validated_domain = new;
Basic Usage
use ;
use SortedSet;
use Deref; // Needed for .deref()
use TryNew; // Needed for .try_new()
// Create uniform grid with equal spacing
let domain = new;
let uniform_grid = uniform;
assert_eq!;
// Create non-uniform grid with custom spacing
let coords = from_unsorted;
let non_uniform_grid = try_from_sorted.unwrap;
// Both grids implement the same IntervalPartition trait
assert_eq!;
assert_eq!;
// Efficient point location
let interval_id = uniform_grid.find_interval_id_of_point;
let interval = uniform_grid.interval;
Working with Different Scalar Types
use ;
use ;
use TryNew;
// Performance-optimal with debug safety (recommended)
type OptimalReal = RealNative64StrictFiniteInDebug;
let optimal_grid = uniform;
// Always validated for safety-critical applications
type SafeReal = RealNative64StrictFinite;
let safe_grid = uniform;
Arbitrary Precision with Rug
For applications requiring arbitrary precision arithmetic:
// Add to Cargo.toml:
// grid1d = { version = "0.1.0", features = ["rug"] }
use ;
use RealRugStrictFinite;
use TryNew;
// 256-bit precision arithmetic
type HighPrecisionReal = ;
let high_precision_grid = uniform;
// All grid operations work seamlessly with arbitrary precision
let hp_1 = try_new.unwrap;
let hp_3 = try_new.unwrap;
let point = hp_1 / hp_3; // 1. / 3.
let interval_id = high_precision_grid.find_interval_id_of_point;
๐ฌ Advanced Features
Adaptive Mesh Refinement
use ;
use BTreeMap; // Needed for refinement plan
use TryNew; // Needed for .try_new()
let base_grid = uniform;
// Uniform refinement: double resolution everywhere
let uniform_refinement = base_grid.refine_uniform;
assert_eq!;
// Selective refinement: refine only specific intervals
let selective_plan = from;
let selective_refinement = uniform_refinement.into_refined_grid.refine;
Grid Union for Multi-Physics
use ;
use SortedSet;
use TryNew; // Needed for .try_new()
let domain = new;
// Physics A: coarse global grid
let grid_a = uniform;
// Physics B: fine local grid
let coords_b = from_unsorted;
let grid_b = try_from_sorted.unwrap;
// Create unified grid preserving mappings to both original grids
let union = try_new.unwrap;
println!;
// Map data between grids
for in union.iter_interval_mappings
Working with Different Interval Types
use *;
let closed = new; // [0, 1]
let open = new; // (0, 1)
let half_open = new; // [0, 1)
let unbounded = new; // [0, +โ)
// All implement IntervalTrait for generic operations
assert!;
assert!;
๐๏ธ 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 ;
use Deref; // Needed for accessing coords
use TryNew; // Needed for .try_new()
let grid = uniform;
let coords = grid.coords.deref;
let dx = coords - coords;
// Initialize solution vector
let mut solution = vec!;
// ... (set initial/boundary conditions)
// Apply finite difference stencils (interior points only)
for i in 1..coords.len-1
Spectral Methods
use ;
use TryNew; // Needed for .try_new()
use Deref; // Needed for accessing coords
// Periodic domain for Fourier spectral methods
let domain = new;
let grid = uniform;
// Perfect for FFT algorithms
let dx = grid.coords - grid.coords;
assert!;
Boundary Layer Refinement
use ;
use SortedSet;
// Create boundary layer mesh with fine spacing near walls
let cfd_grid = create_boundary_layer_mesh;
println!;
๐ ๏ธ Development
Building
# Standard build (uses nightly automatically)
# With arbitrary precision support
# Run tests
# Run tests with arbitrary precision
# Generate documentation
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 theruglibrary- Enable with
--features=rug - Provides
num_valid::RealRugStrictFinite<N>types for N-bit precision
- Enable with
๐งช Testing
The library includes comprehensive test suites:
# Run all tests
# Run tests with arbitrary precision
# Run property-based tests
# Run benchmarks
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
-
Clone the repository:
-
Ensure you have the required dependencies:
- Rust nightly toolchain
-
Run the test suite:
-
Check formatting and linting:
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:
- ๐ Complete Tutorial: Step-by-step guide from basics to advanced usage
- ๐ฏ Decision Guide: Choose the right grid types and scalar types for your needs
- ๐ API Documentation: Complete API reference with examples
- ๐ง Copilot Instructions: Quick reference for AI-assisted development
- ๐ Changelog: Version history and release notes
Quick Links
- Getting Started: See Your First Grid Tutorial
- Which Grid Type?: See Grid Type Decision Tree
- Which Scalar Type?: See Scalar Type Decision Tree
- Real-World Examples: See Applications Section
- Performance Tuning: See Performance Guide
๐ฏ 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.