# Grid1D Decision Guide
This guide helps you make informed choices about which grid types, scalar types, and features to use in your application.
## Quick Decision Trees
### ๐ฏ Which Grid Type Should I Use?
```
START: What kind of spacing do I need?
โ
โโ Equal spacing everywhere?
โ โโ โ
Use Grid1DUniform
โ โข O(1) point location
โ โข Maximum performance
โ โข Perfect for regular discretizations
โ
โโ Custom/adaptive spacing?
โ โโ โ
Use Grid1DNonUniform
โ โข O(log n) point location
โ โข Flexible point distribution
โ โข Ideal for boundary layers, shock capturing
โ
โโ Combining multiple grids?
โโ โ
Use Grid1DUnion
โข Combines grids from different physics
โข Maintains mappings to original grids
โข For multi-physics simulations
```
### ๐ Which Interval Type Should I Use?
```
START: What are your domain boundaries?
โ
โโ Both endpoints are part of the domain?
โ โโ โ
Use IntervalClosed [a, b]
โ โข Most common choice
โ โข Points at a and b are included
โ โข Standard for boundary value problems
โ
โโ Neither endpoint is part of the domain?
โ โโ โ
Use IntervalOpen (a, b)
โ โข Points strictly between a and b
โ โข For open sets, strict inequalities
โ โข Rare in practice
โ
โโ Periodic or partition-friendly boundary?
โ โโ Left closed, right open?
โ โ โโ โ
Use IntervalLowerClosedUpperOpen [a, b)
โ โ โข Common for periodic domains
โ โ โข Non-overlapping partitions
โ โ โข Example: [0, 2ฯ) for angles
โ โ
โ โโ Left open, right closed?
โ โโ โ
Use IntervalLowerOpenUpperClosed (a, b]
โ โข Asymmetric boundaries
โ โข Less common
โ
โโ Semi-infinite domain?
โ โโ Bounded below only?
โ โ โโ โ
Use IntervalLowerClosedUpperUnbounded [a, +โ)
โ โ โข Example: time t โฅ 0
โ โ โข Example: distance r โฅ 0
โ โ
โ โโ Bounded above only?
โ โโ โ
Use IntervalLowerUnboundedUpperClosed (-โ, b]
โ โข Example: x โค 0
โ
โโ Single point (degenerate)?
โโ โ
Use IntervalSingleton {a}
โข Zero-measure domain
โข Special case handling
```
### ๐ข Which Scalar Type Should I Use?
```
START: What are your priorities?
โ
โโ Maximum performance, fully trusted inputs?
โ โโ โ
Use f64
โ โข Zero overhead
โ โข Standard IEEE 754
โ โข โ ๏ธ No validation - can propagate NaN/Inf
โ
โโ Good performance + safety during development?
โ โโ โ
Use RealNative64StrictFiniteInDebug (โญ RECOMMENDED)
โ โข Same as f64 in release builds
โ โข Validates in debug builds
โ โข Perfect development โ production workflow
โ
โโ Safety-critical application?
โ โโ โ
Use RealNative64StrictFinite
โ โข Always validates (small overhead ~5-10%)
โ โข Guaranteed finite values
โ โข For medical, aerospace, financial
โ
โโ Need arbitrary precision?
โโ โ
Use RealRugStrictFinite<N>
โข N-bit precision (e.g., 256, 512, 1024)
โข For scientific computing, symbolic math
โข Requires `--features=rug`
```
## Detailed Comparison Tables
### Grid Type Comparison
| **Point Location** | O(1) analytical | O(log n) binary search | O(log n) on unified grid |
| **Memory Usage** | 8n + 40 bytes | 8n + 40 bytes | 8(n+m) + overhead |
| **Construction** | O(n) | O(n) | O(n+m) |
| **Best For** | Regular meshes, FDM | Adaptive meshes, FEM | Multi-physics coupling |
| **Spacing** | Uniform | Arbitrary | Combined |
| **Flexibility** | Low | High | Very High |
| **Performance** | โกโกโก Excellent | โกโก Good | โกโก Good |
*Where n = number of points, m = points in second grid*
### Scalar Type Comparison
| **f64** | โ None | โกโกโก | โกโกโก | Trusted inputs, maximum speed |
| **RealNative64StrictFiniteInDebug** | โ
Debug only | โกโก | โกโกโก | **Recommended default** |
| **RealNative64StrictFinite** | โ
Always | โกโก | โกโก | Safety-critical |
| **RealRugStrictFinite\<N\>** | โ
Always | โก | โก | Arbitrary precision |
### Interval Type Selection
| Include both endpoints | `IntervalClosed` | [a, b] | Standard bounded domains |
| Exclude both endpoints | `IntervalOpen` | (a, b) | Open sets, strict inequalities |
| Include left, exclude right | `IntervalLowerClosedUpperOpen` | [a, b) | Periodic domains, half-open ranges |
| Include right, exclude left | `IntervalLowerOpenUpperClosed` | (a, b] | Asymmetric boundaries |
| Unbounded above | `IntervalLowerClosedUpperUnbounded` | [a, +โ) | Semi-infinite domains |
| Unbounded below | `IntervalLowerUnboundedUpperClosed` | (-โ, b] | Semi-infinite domains |
| Single point | `IntervalSingleton` | {a} | Degenerate intervals |
### ๐ Which Refinement Strategy Should I Use?
```
START: How do you want to refine your grid?
โ
โโ Same refinement everywhere?
โ โโ โ
Use refine_uniform()
โ โข Every interval gets the same number of extra points
โ โข Simple API: grid.refine_uniform(&num_extra_points)
โ โข Best for: uniform quality improvement
โ
โโ Different refinement in different regions?
โ โโ โ
Use refine() with BTreeMap
โ โข Specify extra points per interval ID
โ โข Only refine where needed
โ โข Best for: adaptive mesh refinement (AMR)
โ
โโ Combining two existing grids?
โโ โ
Use Grid1DUnion instead
โข Not refinement, but unification
โข Preserves both original grid structures
โข Best for: multi-physics coupling
```
### ๐ Grid Union vs Grid Refinement: When to Use Each?
```
START: What is your goal?
โ
โโ Make current grid finer?
โ โ
โ โโ Uniformly finer everywhere?
โ โ โโ โ
Use refine_uniform()
โ โ
โ โโ Finer only in specific regions?
โ โโ โ
Use refine() with BTreeMap
โ
โโ Couple two independent grids?
โ โโ โ
Use Grid1DUnion
โ โข Both grids retain their identity
โ โข Bidirectional mapping preserved
โ โข For multi-physics data exchange
โ
โโ Create grid from scratch with custom points?
โโ โ
Use Grid1D::try_from_sorted()
โข Direct construction from coordinates
โข Maximum flexibility
โข When refinement hierarchy not needed
```
| **Purpose** | Global refinement | Selective refinement | Grid coupling |
| **Tracks parent intervals** | โ
Yes | โ
Yes | โ
Yes (both grids) |
| **Preserves original grid** | As parent | As parent | Both grids intact |
| **Best for** | h-refinement FEM | AMR, error-driven | Multi-physics |
| **Memory overhead** | Low | Low | Medium |
| **Complexity** | Simple | Medium | Higher |
### Refinement Strategy Examples
**Uniform refinement** - Every interval gets 2 extra points (โ 3 sub-intervals each):
```rust
use grid1d::{Grid1D, IntervalPartition, intervals::*, scalars::{NumIntervals, PositiveNumPoints1D}};
use try_create::TryNew;
let grid = Grid1D::uniform(IntervalClosed::new(0.0, 1.0), NumIntervals::try_new(4).unwrap());
// Original: 4 intervals
let refinement = grid.refine_uniform(&PositiveNumPoints1D::try_new(2).unwrap());
// Refined: 12 intervals (4 ร 3)
// Track parent-child relationships
for (refined_id, original_id) in refinement.iter_refined_with_mapping() {
println!("Refined interval {} came from original interval {}",
refined_id.as_ref(), original_id.as_ref());
}
```
**Selective refinement** - Only refine where needed (adaptive mesh refinement):
```rust
use grid1d::{Grid1D, IntervalPartition, intervals::*, scalars::{NumIntervals, PositiveNumPoints1D, IntervalId}};
use std::collections::BTreeMap;
use try_create::TryNew;
let grid = Grid1D::uniform(IntervalClosed::new(0.0, 1.0), NumIntervals::try_new(4).unwrap());
// Refine only intervals 0 and 2
let plan = BTreeMap::from([
(IntervalId::new(0), PositiveNumPoints1D::try_new(3).unwrap()), // 4 sub-intervals
(IntervalId::new(2), PositiveNumPoints1D::try_new(1).unwrap()), // 2 sub-intervals
]);
let refinement = grid.refine(&plan);
// Result: intervals 1 and 3 unchanged, 0 and 2 refined
```
**Grid union** - Couple two grids for multi-physics:
```rust
use grid1d::{Grid1D, Grid1DUnion, IntervalPartition, intervals::*, scalars::NumIntervals};
use try_create::TryNew;
// Flow solver: coarse grid
let flow_grid = Grid1D::uniform(
IntervalClosed::new(0.0, 1.0),
NumIntervals::try_new(10).unwrap()
);
// Chemistry solver: fine grid in reaction zone
let chem_grid = Grid1D::uniform(
IntervalClosed::new(0.0, 1.0),
NumIntervals::try_new(50).unwrap()
);
// Unified grid for coupling
let coupled = Grid1DUnion::try_new(&flow_grid, &chem_grid).unwrap();
// Transfer data between physics
for (unified_id, flow_id, chem_id) in coupled.iter_interval_mappings() {
// Map solution from flow grid to chemistry grid and vice versa
}
```
## Use Case โ Recommended Configuration
### Finite Difference Methods (Regular Mesh)
**Scenario**: Solving PDEs on uniform grids
```rust
// Recommended configuration
type Real = RealNative64StrictFiniteInDebug; // Debug safety
type Domain = IntervalClosed<Real>; // Closed boundaries
let grid = Grid1D::uniform(
Domain::new(Real::try_new(0.0).unwrap(), Real::try_new(1.0).unwrap()),
NumIntervals::try_new(1000).unwrap()
);
```
**Why?**
- โ
Uniform grid: O(1) point location for maximum performance
- โ
Debug validation: Catches NaN/Inf during development
- โ
Closed interval: Standard for boundary value problems
- โ
Release performance: Same as raw f64
### Finite Element Methods (Adaptive Mesh)
**Scenario**: FEM with mesh refinement in regions of interest
```rust
use grid1d::{Grid1D, IntervalPartition, intervals::*, scalars::{NumIntervals, PositiveNumPoints1D, IntervalId}};
use num_valid::RealNative64StrictFinite;
use std::collections::BTreeMap;
use try_create::TryNew;
// Recommended configuration
type Real = RealNative64StrictFinite; // Always validated
type Domain = IntervalClosed<Real>;
// Start with coarse base mesh
let base_grid = Grid1D::uniform(
Domain::new(Real::try_new(0.0).unwrap(), Real::try_new(1.0).unwrap()),
NumIntervals::try_new(10).unwrap()
);
// Compute error estimator, identify intervals needing refinement
let error_threshold = 0.01;
let refinement_plan: BTreeMap<IntervalId, PositiveNumPoints1D> = (0..*base_grid.num_intervals().as_ref())
.filter(|&i| estimate_error(i) > error_threshold)
.map(|i| (IntervalId::new(i), PositiveNumPoints1D::try_new(1).unwrap()))
.collect();
// Refine selectively where solution changes rapidly
let refined = base_grid.refine(&refinement_plan);
fn estimate_error(_interval: usize) -> f64 { 0.02 } // Placeholder
```
**Why?**
- โ
Selective refinement: Only refine where errors are large
- โ
Always validated: Prevents propagation of numerical errors
- โ
Refinement tracking: Maintains parent-child relationships
- โ
Flexible spacing: Optimal resolution where needed
### Computational Fluid Dynamics (Boundary Layers)
**Scenario**: Viscous flow with wall boundary layers
```rust
use grid1d::{Grid1D, intervals::*};
use num_valid::RealNative64StrictFiniteInDebug;
use sorted_vec::partial::SortedSet;
// Recommended configuration
type Real = RealNative64StrictFiniteInDebug;
type Domain = IntervalClosed<Real>;
// Generate boundary layer mesh with geometric stretching
fn create_boundary_layer_distribution(
wall_position: f64,
outer_position: f64,
first_cell_height: f64,
stretching_ratio: f64,
max_points: usize,
) -> Vec<f64> {
let mut coords = vec![wall_position];
let mut y = wall_position;
let mut dy = first_cell_height;
for _ in 1..max_points {
y += dy;
if y >= outer_position {
coords.push(outer_position);
break;
}
coords.push(y);
dy *= stretching_ratio;
}
coords
}
let coords = create_boundary_layer_distribution(0.0, 1.0, 0.001, 1.2, 50);
// Domain is inferred from first and last coordinates
let grid = Grid1D::<IntervalClosed<f64>>::try_from_sorted(
SortedSet::from_unsorted(coords)
).unwrap();
```
**Why?**
- โ
Non-uniform grid: Fine near wall, coarse in freestream
- โ
Custom distribution: Precise control over point clustering
- โ
Performance: Fast lookups even with non-uniform spacing
- โ
Simplified API: Domain inferred from coordinates
### Multi-Physics Coupling
**Scenario**: Coupling fluid flow with chemical reactions
```rust
use grid1d::{Grid1D, Grid1DUnion, IntervalPartition, intervals::*, scalars::NumIntervals};
use num_valid::RealNative64StrictFinite;
use sorted_vec::partial::SortedSet;
use try_create::TryNew;
// Recommended configuration
type Real = RealNative64StrictFinite; // Safety for complex coupling
type Domain = IntervalClosed<Real>;
// Flow grid: coarse global mesh
let flow_grid = Grid1D::uniform(
Domain::new(Real::try_new(0.0).unwrap(), Real::try_new(1.0).unwrap()),
NumIntervals::try_new(20).unwrap()
);
// Chemistry grid: fine mesh in reaction zone (non-uniform)
let chem_coords: Vec<f64> = (0..=100).map(|i| {
let t = i as f64 / 100.0;
// Cluster points near x=0.5 (reaction zone)
0.5 + 0.5 * (std::f64::consts::PI * (t - 0.5)).tanh()
}).collect();
let chemistry_grid = Grid1D::<IntervalClosed<f64>>::try_from_sorted(
SortedSet::from_unsorted(chem_coords)
).unwrap();
// Unified grid for data transfer
let coupled = Grid1DUnion::try_new(&flow_grid, &chemistry_grid).unwrap();
// Map between grids for data exchange
for (unified_id, flow_id, chem_id) in coupled.iter_interval_mappings() {
// Transfer temperature from flow โ chemistry
// Transfer reaction rates from chemistry โ flow
println!("Unified {} maps to flow {}, chemistry {}",
unified_id.as_ref(), flow_id.as_ref(), chem_id.as_ref());
}
```
**Why?**
- โ
Grid union: Preserves both discretizations
- โ
Bidirectional mapping: Efficient data transfer
- โ
Validated scalars: Prevents coupling instabilities
- โ
Mixed grid types: Uniform flow + non-uniform chemistry
### High-Precision Scientific Computing
**Scenario**: Numerical analysis requiring exact arithmetic
```rust
#[cfg(feature = "rug")]
{
// Recommended configuration
type Real = RealRugStrictFinite<256>; // 256-bit precision
type Domain = IntervalClosed<Real>;
let grid = Grid1D::uniform(
Domain::new(
Real::try_new(0.0).unwrap(),
Real::try_new(1.0).unwrap()
),
NumIntervals::try_new(100).unwrap()
);
}
```
**Why?**
- โ
Arbitrary precision: Eliminates rounding errors
- โ
Validated: Maintains numerical integrity
- โ
Same API: Easy to switch from f64
### Real-Time Applications
**Scenario**: Simulation loop with strict timing requirements
```rust
// Recommended configuration
type Real = f64; // Maximum performance, no validation
type Domain = IntervalClosed<Real>;
let grid = Grid1D::uniform(
Domain::new(0.0, 1.0),
NumIntervals::try_new(1000).unwrap()
);
// In real-time loop
loop {
// O(1) point location - predictable timing
let interval_id = grid.find_interval_id_of_point(&sensor_value);
// ... process ...
}
```
**Why?**
- โ
Raw f64: Zero overhead, predictable performance
- โ
Uniform grid: O(1) lookups, no branch prediction issues
- โ
No validation: No runtime checks in critical path
## Performance Tuning Guidelines
### When to Choose Uniform Grid
โ
**Choose uniform grid when:**
- All intervals need equal spacing
- Maximum performance is critical
- Simple, regular domain discretization
- Point location happens frequently
โ **Avoid uniform grid when:**
- Features occur at different scales
- Adaptive refinement is needed
- Boundary layer resolution required
### When to Choose Non-Uniform Grid
โ
**Choose non-uniform grid when:**
- Adaptive spacing is beneficial
- Multiple length scales present
- Boundary layer capture needed
- Features localized in space
โ **Avoid non-uniform grid when:**
- Uniform spacing is sufficient
- Maximum performance is critical
- Simple regular problem
### When to Use Parallel Point Location
โ
**Use `find_intervals_for_points_parallel()` when:**
- Locating โฅ1,000 points at once
- Running on multi-core systems
- Point location is a bottleneck
- Processing large datasets
โ **Use sequential `find_intervals_for_points()` when:**
- Locating < 1,000 points
- Single-threaded environment required
- Need deterministic ordering
- Rayon overhead would dominate
**Performance guidance:**
| 100 | โ
Faster | โ Overhead | N/A |
| 1,000 | โ Equal | โ Equal | ~1ร |
| 10,000 | Slower | โ
Faster | ~3-4ร |
| 100,000+ | Much slower | โ
Much faster | ~6-8ร |
### When to Use Grid Union
โ
**Use grid union when:**
- Coupling multiple physics
- Different solvers on different grids
- Need to preserve original grid structure
- Frequent data transfer between grids
โ **Avoid grid union when:**
- Single physics problem
- Don't need grid mapping
- Memory is constrained
## Scalar Type Performance Impact
### Benchmarked Overheads (relative to f64)
| Arithmetic | 1.0ร | 1.5-2.0ร | 1.0ร | 1.05-1.10ร |
| Comparison | 1.0ร | 1.3-1.5ร | 1.0ร | 1.02-1.05ร |
| Grid creation | 1.0ร | 1.5ร | 1.0ร | 1.10ร |
| Point location | 1.0ร | 1.2ร | 1.0ร | 1.05ร |
**Key takeaway**: `RealNative64StrictFiniteInDebug` has zero overhead in release builds.
## Migration Strategies
### From f64 to Validated Types
**Step 1**: Start with raw f64
```rust
let grid = Grid1D::uniform(
IntervalClosed::new(0.0_f64, 1.0_f64),
NumIntervals::try_new(100).unwrap()
);
```
**Step 2**: Add debug validation (no code changes needed in release!)
```rust
type Real = RealNative64StrictFiniteInDebug;
let grid = Grid1D::uniform(
IntervalClosed::new(
Real::try_new(0.0).unwrap(),
Real::try_new(1.0).unwrap()
),
NumIntervals::try_new(100).unwrap()
);
```
**Step 3**: If safety-critical, use always-validated
```rust
type Real = RealNative64StrictFinite;
// Same code as Step 2
```
### From Uniform to Non-Uniform
**Before** (uniform):
```rust
use grid1d::{Grid1D, intervals::*, scalars::NumIntervals};
use try_create::TryNew;
let grid = Grid1D::uniform(
IntervalClosed::new(0.0, 1.0),
NumIntervals::try_new(100).unwrap()
);
```
**After** (non-uniform with custom spacing):
```rust
use grid1d::{Grid1D, intervals::*};
use sorted_vec::partial::SortedSet;
fn generate_custom_distribution() -> Vec<f64> {
// Example: cluster points near boundaries
(0..=100).map(|i| {
let t = i as f64 / 100.0;
// Chebyshev-like distribution
0.5 * (1.0 - (std::f64::consts::PI * t).cos())
}).collect()
}
let coords = generate_custom_distribution();
// Domain is inferred from first and last coordinates
let grid = Grid1D::<IntervalClosed<f64>>::try_from_sorted(
SortedSet::from_unsorted(coords)
).unwrap();
```
## Common Pitfalls and Solutions
### Pitfall 1: Using Wrong Grid Type
โ **Wrong**:
```rust
use grid1d::{Grid1D, intervals::*};
use sorted_vec::partial::SortedSet;
// Using non-uniform grid when uniform would work
SortedSet::from_unsorted(coords)
).unwrap(); // O(log n) lookups - unnecessary overhead!
```
โ
**Correct**:
```rust
use grid1d::{Grid1D, intervals::*, scalars::NumIntervals};
use try_create::TryNew;
// Use uniform grid for better performance
let grid = Grid1D::uniform(
IntervalClosed::new(0.0, 1.0),
NumIntervals::try_new(100).unwrap()
); // O(1) lookups!
```
### Pitfall 2: Over-validating in Performance-Critical Code
โ **Wrong**:
```rust
type Real = RealNative64StrictFinite; // Always validates
// In hot loop
for _ in 0..1_000_000 {
let result = a + b; // Validates every operation
}
```
โ
**Correct**:
```rust
type Real = RealNative64StrictFiniteInDebug; // Validates only in debug
// In hot loop - no overhead in release
for _ in 0..1_000_000 {
let result = a + b; // Fast in release
}
```
### Pitfall 3: Forgetting to Handle Errors
โ **Wrong**:
```rust
use grid1d::{Grid1D, intervals::*};
use sorted_vec::partial::SortedSet;
let coords = SortedSet::from_unsorted(vec![0.0, 0.5, 1.0]);
let grid = Grid1D::<IntervalClosed<f64>>::try_from_sorted(coords).unwrap(); // Panics on error!
```
โ
**Correct**:
```rust
use grid1d::{Grid1D, intervals::*};
use sorted_vec::partial::SortedSet;
let coords = SortedSet::from_unsorted(vec![0.0, 0.5, 1.0]);
let grid = Grid1D::<IntervalClosed<f64>>::try_from_sorted(coords)
.map_err(|e| {
eprintln!("Grid creation failed: {:?}", e);
e
})?;
```
### Pitfall 4: Using new() with Untrusted Input
โ **Wrong**:
```rust
use grid1d::intervals::*;
// new() panics if lower > upper!
let user_input_lo = get_user_input(); // Could be anything
let user_input_hi = get_user_input();
let interval = IntervalClosed::new(user_input_lo, user_input_hi); // ๐ฅ May panic!
fn get_user_input() -> f64 { 5.0 } // Placeholder
```
โ
**Correct**:
```rust
use grid1d::intervals::*;
let user_input_lo = get_user_input();
let user_input_hi = get_user_input();
let interval = IntervalClosed::try_new(user_input_lo, user_input_hi)
.map_err(|e| format!("Invalid interval: {:?}", e))?;
fn get_user_input() -> f64 { 5.0 } // Placeholder
```
## Summary Recommendations
### For Most Users (โญ Recommended)
```rust
use grid1d::{Grid1D, intervals::*, scalars::NumIntervals};
use num_valid::RealNative64StrictFiniteInDebug;
use try_create::TryNew;
type Real = RealNative64StrictFiniteInDebug;
type Domain = IntervalClosed<Real>;
// Uniform grid for regular problems
let grid = Grid1D::uniform(
Domain::new(Real::try_new(0.0).unwrap(), Real::try_new(1.0).unwrap()),
NumIntervals::try_new(100).unwrap()
);
```
**This gives you**:
- โ
Maximum performance in release builds
- โ
Safety validation during development
- โ
O(1) point location for uniform grids
- โ
Easy to change if needs evolve
### Quick Reference Card
| Regular mesh, maximum speed | `Grid1DUniform` + `f64` |
| **Most applications** | **`Grid1DUniform` + `RealNative64StrictFiniteInDebug`** |
| Adaptive mesh | `Grid1DNonUniform` + `RealNative64StrictFiniteInDebug` |
| Safety-critical | Any grid + `RealNative64StrictFinite` |
| Multi-physics | `Grid1DUnion` + `RealNative64StrictFinite` |
| High precision | Any grid + `RealRugStrictFinite<N>` |
| Periodic domains | Any grid + `IntervalLowerClosedUpperOpen` |
---
For more details, see:
- [Complete Tutorial](./TUTORIAL.md)
- [API Documentation](https://docs.rs/grid1d)
- [README](./README.md)