Skip to main content

ruvector_solver/
neumann.rs

1//! Jacobi-preconditioned Neumann Series iterative solver.
2//!
3//! Solves the linear system `Ax = b` by splitting `A = D - R` (where `D` is
4//! the diagonal part) and iterating:
5//!
6//! ```text
7//! x_{k+1} = x_k + D^{-1} (b - A x_k)
8//! ```
9//!
10//! This is equivalent to the Neumann series `x = sum_{k=0}^{K} M^k D^{-1} b`
11//! where `M = I - D^{-1} A`. Convergence requires `rho(M) < 1`, which is
12//! guaranteed for strictly diagonally dominant matrices.
13//!
14//! # Algorithm
15//!
16//! The iteration maintains a running solution `x` and residual `r = b - Ax`:
17//!
18//! ```text
19//! x_0 = D^{-1} b
20//! for k = 0, 1, 2, ...:
21//!     r = b - A * x_k
22//!     x_{k+1} = x_k + D^{-1} * r
23//!     if ||r|| < tolerance:
24//!         break
25//! ```
26//!
27//! # Convergence
28//!
29//! Before solving, the solver estimates `rho(I - D^{-1}A)` via a 10-step
30//! power iteration and rejects the problem with
31//! [`SolverError::SpectralRadiusExceeded`] if `rho >= 1.0`. During iteration,
32//! if the residual grows by more than 2x between consecutive steps,
33//! [`SolverError::NumericalInstability`] is returned.
34
35use std::time::Instant;
36
37use tracing::{debug, info, instrument, warn};
38
39use crate::error::{SolverError, ValidationError};
40use crate::traits::SolverEngine;
41use crate::types::{
42    Algorithm, ComplexityClass, ComplexityEstimate, ComputeBudget, ConvergenceInfo, CsrMatrix,
43    SolverResult, SparsityProfile,
44};
45
46/// Number of power-iteration steps used to estimate the spectral radius.
47const POWER_ITERATION_STEPS: usize = 10;
48
49/// If the residual grows by more than this factor in a single step, the solver
50/// declares numerical instability.
51const INSTABILITY_GROWTH_FACTOR: f64 = 2.0;
52
53// ---------------------------------------------------------------------------
54// NeumannSolver
55// ---------------------------------------------------------------------------
56
57/// Neumann Series solver for sparse linear systems.
58///
59/// Computes `x = sum_{k=0}^{K} (I - A)^k * b` by maintaining a residual
60/// vector and accumulating partial sums until convergence.
61///
62/// # Example
63///
64/// ```rust
65/// use ruvector_solver::types::CsrMatrix;
66/// use ruvector_solver::neumann::NeumannSolver;
67///
68/// // Diagonally dominant 2x2: A = [[2, -0.5], [-0.5, 2]]
69/// let a = CsrMatrix::<f32>::from_coo(2, 2, vec![
70///     (0, 0, 2.0_f32), (0, 1, -0.5_f32),
71///     (1, 0, -0.5_f32), (1, 1, 2.0_f32),
72/// ]);
73/// let b = vec![1.0_f32, 1.0];
74///
75/// let solver = NeumannSolver::new(1e-6, 500);
76/// let result = solver.solve(&a, &b).unwrap();
77/// assert!(result.residual_norm < 1e-4);
78/// ```
79#[derive(Debug, Clone)]
80pub struct NeumannSolver {
81    /// Target residual L2 norm for convergence.
82    pub tolerance: f64,
83    /// Maximum number of iterations before giving up.
84    pub max_iterations: usize,
85}
86
87impl NeumannSolver {
88    /// Create a new `NeumannSolver`.
89    ///
90    /// # Arguments
91    ///
92    /// * `tolerance` - Stop when `||r|| < tolerance`.
93    /// * `max_iterations` - Upper bound on iterations.
94    pub fn new(tolerance: f64, max_iterations: usize) -> Self {
95        Self {
96            tolerance,
97            max_iterations,
98        }
99    }
100
101    /// Estimate the spectral radius of `M = I - D^{-1}A` via 10-step power
102    /// iteration.
103    ///
104    /// Runs [`POWER_ITERATION_STEPS`] iterations of the power method on the
105    /// Jacobi iteration matrix `M = I - D^{-1}A`. Returns the Rayleigh-quotient
106    /// estimate of the dominant eigenvalue magnitude.
107    ///
108    /// # Arguments
109    ///
110    /// * `matrix` - The coefficient matrix `A` (must be square).
111    ///
112    /// # Returns
113    ///
114    /// Estimated `|lambda_max(I - D^{-1}A)|`. If this is `>= 1.0`, the
115    /// Jacobi-preconditioned Neumann series will diverge.
116    #[instrument(skip(matrix), fields(n = matrix.rows))]
117    pub fn estimate_spectral_radius(matrix: &CsrMatrix<f32>) -> f64 {
118        let n = matrix.rows;
119        if n == 0 {
120            return 0.0;
121        }
122
123        let d_inv = extract_diag_inv_f32(matrix);
124        Self::estimate_spectral_radius_with_diag(matrix, &d_inv)
125    }
126
127    /// Inner helper: estimate spectral radius using a pre-computed `d_inv`.
128    ///
129    /// This avoids recomputing the diagonal inverse when the caller already
130    /// has it (e.g. `solve()` needs `d_inv` for both the spectral check and
131    /// the Jacobi iteration).
132    fn estimate_spectral_radius_with_diag(matrix: &CsrMatrix<f32>, d_inv: &[f32]) -> f64 {
133        let n = matrix.rows;
134        if n == 0 {
135            return 0.0;
136        }
137
138        // Initialise with a deterministic pseudo-random unit vector.
139        let mut v: Vec<f32> = (0..n)
140            .map(|i| ((i.wrapping_mul(7).wrapping_add(13)) % 100) as f32 / 100.0)
141            .collect();
142        let norm = l2_norm_f32(&v);
143        if norm > 1e-12 {
144            scale_vec_f32(&mut v, 1.0 / norm);
145        }
146
147        let mut av = vec![0.0f32; n]; // scratch for A*v
148        let mut w = vec![0.0f32; n]; // scratch for M*v = v - D^{-1}*A*v
149        let mut eigenvalue_estimate = 0.0f64;
150
151        for _ in 0..POWER_ITERATION_STEPS {
152            // w = v - D^{-1} * A * v  (i.e. M * v)
153            matrix.spmv(&v, &mut av);
154            for j in 0..n {
155                w[j] = v[j] - d_inv[j] * av[j];
156            }
157
158            // Rayleigh quotient: lambda = v^T w  (v is unit-length).
159            let dot: f64 = v
160                .iter()
161                .zip(w.iter())
162                .map(|(&a, &b)| a as f64 * b as f64)
163                .sum();
164            eigenvalue_estimate = dot;
165
166            // Normalise w -> v for the next step.
167            let w_norm = l2_norm_f32(&w);
168            if w_norm < 1e-12 {
169                break;
170            }
171            for j in 0..n {
172                v[j] = w[j] / w_norm as f32;
173            }
174        }
175
176        let rho = eigenvalue_estimate.abs();
177        debug!(rho, "estimated spectral radius of (I - D^-1 A)");
178        rho
179    }
180
181    /// Core Jacobi-preconditioned Neumann-series solve operating on `f32`.
182    ///
183    /// Validates inputs, checks the spectral radius of `I - D^{-1}A` via
184    /// power iteration, then runs the iteration returning a [`SolverResult`].
185    ///
186    /// # Errors
187    ///
188    /// - [`SolverError::InvalidInput`] if the matrix is non-square or the RHS
189    ///   length does not match.
190    /// - [`SolverError::SpectralRadiusExceeded`] if `rho(I - D^{-1}A) >= 1`.
191    /// - [`SolverError::NumericalInstability`] if the residual grows by more
192    ///   than 2x in a single step.
193    /// - [`SolverError::NonConvergence`] if the iteration budget is exhausted.
194    #[instrument(skip(self, matrix, rhs), fields(n = matrix.rows, nnz = matrix.nnz()))]
195    pub fn solve(&self, matrix: &CsrMatrix<f32>, rhs: &[f32]) -> Result<SolverResult, SolverError> {
196        let start = Instant::now();
197        let n = matrix.rows;
198
199        // ------------------------------------------------------------------
200        // Input validation
201        // ------------------------------------------------------------------
202        if matrix.rows != matrix.cols {
203            return Err(SolverError::InvalidInput(
204                ValidationError::DimensionMismatch(format!(
205                    "matrix must be square: got {}x{}",
206                    matrix.rows, matrix.cols,
207                )),
208            ));
209        }
210
211        if rhs.len() != n {
212            return Err(SolverError::InvalidInput(
213                ValidationError::DimensionMismatch(format!(
214                    "rhs length {} does not match matrix dimension {}",
215                    rhs.len(),
216                    n,
217                )),
218            ));
219        }
220
221        // Edge case: empty system.
222        if n == 0 {
223            return Ok(SolverResult {
224                solution: Vec::new(),
225                iterations: 0,
226                residual_norm: 0.0,
227                wall_time: start.elapsed(),
228                convergence_history: Vec::new(),
229                algorithm: Algorithm::Neumann,
230            });
231        }
232
233        // Extract D^{-1} once — reused for both the spectral radius check
234        // and the Jacobi-preconditioned iteration that follows.
235        let d_inv = extract_diag_inv_f32(matrix);
236
237        // ------------------------------------------------------------------
238        // Spectral radius pre-check (10-step power iteration on I - D^{-1}A)
239        // ------------------------------------------------------------------
240        let rho = Self::estimate_spectral_radius_with_diag(matrix, &d_inv);
241        if rho >= 1.0 {
242            warn!(rho, "spectral radius >= 1.0, Neumann series will diverge");
243            return Err(SolverError::SpectralRadiusExceeded {
244                spectral_radius: rho,
245                limit: 1.0,
246                algorithm: Algorithm::Neumann,
247            });
248        }
249        info!(rho, "spectral radius check passed");
250
251        // ------------------------------------------------------------------
252        // Jacobi-preconditioned iteration (fused kernel)
253        //
254        //   x_0 = D^{-1} * b
255        //   loop:
256        //       r = b - A * x_k            (fused with norm computation)
257        //       if ||r|| < tolerance: break
258        //       x_{k+1} = x_k + D^{-1} * r (fused with residual storage)
259        //
260        // Key optimization: uses fused_residual_norm_sq to compute
261        // r = b - Ax and ||r||^2 in a single pass, avoiding a separate
262        // spmv + subtraction + norm computation (3 memory traversals -> 1).
263        // ------------------------------------------------------------------
264        let mut x: Vec<f32> = (0..n).map(|i| d_inv[i] * rhs[i]).collect();
265        let mut r = vec![0.0f32; n]; // residual buffer (reused each iteration)
266
267        let mut convergence_history = Vec::with_capacity(self.max_iterations.min(256));
268        let mut prev_residual_norm = f64::MAX;
269        let final_residual_norm: f64;
270        let mut iterations_done: usize = 0;
271
272        for k in 0..self.max_iterations {
273            // Fused: compute r = b - Ax and ||r||^2 in one pass.
274            let residual_norm_sq = matrix.fused_residual_norm_sq(&x, rhs, &mut r);
275            let residual_norm = residual_norm_sq.sqrt();
276            iterations_done = k + 1;
277
278            convergence_history.push(ConvergenceInfo {
279                iteration: k,
280                residual_norm,
281            });
282
283            debug!(iteration = k, residual_norm, "neumann iteration");
284
285            // Convergence check.
286            if residual_norm < self.tolerance {
287                final_residual_norm = residual_norm;
288                info!(iterations = iterations_done, residual_norm, "converged");
289                return Ok(SolverResult {
290                    solution: x,
291                    iterations: iterations_done,
292                    residual_norm: final_residual_norm,
293                    wall_time: start.elapsed(),
294                    convergence_history,
295                    algorithm: Algorithm::Neumann,
296                });
297            }
298
299            // NaN / Inf guard.
300            if residual_norm.is_nan() || residual_norm.is_infinite() {
301                return Err(SolverError::NumericalInstability {
302                    iteration: k,
303                    detail: format!("residual became {residual_norm}"),
304                });
305            }
306
307            // Instability check: residual grew by > 2x.
308            if k > 0
309                && prev_residual_norm < f64::MAX
310                && prev_residual_norm > 0.0
311                && residual_norm > INSTABILITY_GROWTH_FACTOR * prev_residual_norm
312            {
313                warn!(
314                    iteration = k,
315                    prev = prev_residual_norm,
316                    current = residual_norm,
317                    "residual diverging",
318                );
319                return Err(SolverError::NumericalInstability {
320                    iteration: k,
321                    detail: format!(
322                        "residual grew from {prev_residual_norm:.6e} to \
323                         {residual_norm:.6e} (>{INSTABILITY_GROWTH_FACTOR:.0}x)",
324                    ),
325                });
326            }
327
328            // Fused update: x[j] += d_inv[j] * r[j]
329            // 4-wide unrolled for ILP.
330            let chunks = n / 4;
331            for c in 0..chunks {
332                let j = c * 4;
333                x[j] += d_inv[j] * r[j];
334                x[j + 1] += d_inv[j + 1] * r[j + 1];
335                x[j + 2] += d_inv[j + 2] * r[j + 2];
336                x[j + 3] += d_inv[j + 3] * r[j + 3];
337            }
338            for j in (chunks * 4)..n {
339                x[j] += d_inv[j] * r[j];
340            }
341
342            prev_residual_norm = residual_norm;
343        }
344
345        // Exhausted iteration budget without converging.
346        final_residual_norm = prev_residual_norm;
347        Err(SolverError::NonConvergence {
348            iterations: iterations_done,
349            residual: final_residual_norm,
350            tolerance: self.tolerance,
351        })
352    }
353}
354
355// ---------------------------------------------------------------------------
356// SolverEngine trait implementation (f64 interface)
357// ---------------------------------------------------------------------------
358
359impl SolverEngine for NeumannSolver {
360    /// Solve via the Neumann series.
361    ///
362    /// Adapts the `f64` trait interface to the internal `f32` solver by
363    /// converting the input matrix and RHS, running the solver, then
364    /// returning the `f32` solution.
365    fn solve(
366        &self,
367        matrix: &CsrMatrix<f64>,
368        rhs: &[f64],
369        budget: &ComputeBudget,
370    ) -> Result<SolverResult, SolverError> {
371        let start = Instant::now();
372
373        // Validate that f64 values fit in f32 range.
374        for (i, &v) in matrix.values.iter().enumerate() {
375            if v.is_finite() && v.abs() > f32::MAX as f64 {
376                return Err(SolverError::InvalidInput(ValidationError::NonFiniteValue(
377                    format!("matrix value at index {i} ({v:.6e}) overflows f32"),
378                )));
379            }
380        }
381        for (i, &v) in rhs.iter().enumerate() {
382            if v.is_finite() && v.abs() > f32::MAX as f64 {
383                return Err(SolverError::InvalidInput(ValidationError::NonFiniteValue(
384                    format!("rhs value at index {i} ({v:.6e}) overflows f32"),
385                )));
386            }
387        }
388
389        // Convert f64 matrix to f32 for the core solver.
390        let f32_matrix = CsrMatrix {
391            row_ptr: matrix.row_ptr.clone(),
392            col_indices: matrix.col_indices.clone(),
393            values: matrix.values.iter().map(|&v| v as f32).collect(),
394            rows: matrix.rows,
395            cols: matrix.cols,
396        };
397        let f32_rhs: Vec<f32> = rhs.iter().map(|&v| v as f32).collect();
398
399        // Use the tighter of the solver's own tolerance and the caller's budget,
400        // but no tighter than f32 precision allows (the Neumann solver operates
401        // internally in f32, so residuals below ~f32::EPSILON are unreachable).
402        let max_iters = self.max_iterations.min(budget.max_iterations);
403        let tol = self
404            .tolerance
405            .min(budget.tolerance)
406            .max(f32::EPSILON as f64 * 4.0);
407
408        let inner_solver = NeumannSolver::new(tol, max_iters);
409
410        let mut result = inner_solver.solve(&f32_matrix, &f32_rhs)?;
411
412        // Check wall-time budget.
413        if start.elapsed() > budget.max_time {
414            return Err(SolverError::BudgetExhausted {
415                reason: "wall-clock time limit exceeded".to_string(),
416                elapsed: start.elapsed(),
417            });
418        }
419
420        // Adjust wall time to include conversion overhead.
421        result.wall_time = start.elapsed();
422        Ok(result)
423    }
424
425    fn estimate_complexity(&self, profile: &SparsityProfile, n: usize) -> ComplexityEstimate {
426        // Estimated iterations: ceil( ln(1/tol) / |ln(rho)| )
427        let rho = profile.estimated_spectral_radius.max(0.01).min(0.999);
428        let est_iters = ((1.0 / self.tolerance).ln() / (1.0 - rho).ln().abs()).ceil() as usize;
429        let est_iters = est_iters.min(self.max_iterations).max(1);
430
431        ComplexityEstimate {
432            algorithm: Algorithm::Neumann,
433            // Each iteration does one SpMV (2 * nnz flops) + O(n) vector ops.
434            estimated_flops: (est_iters as u64) * (profile.nnz as u64) * 2,
435            estimated_iterations: est_iters,
436            // Working memory: x, r, ar (3 vectors of f32).
437            estimated_memory_bytes: n * 4 * 3,
438            complexity_class: ComplexityClass::SublinearNnz,
439        }
440    }
441
442    fn algorithm(&self) -> Algorithm {
443        Algorithm::Neumann
444    }
445}
446
447// ---------------------------------------------------------------------------
448// Internal helpers
449// ---------------------------------------------------------------------------
450
451/// Extract `D^{-1}` from a CSR matrix (the reciprocal of each diagonal entry).
452///
453/// If a diagonal entry is zero or very small, uses `1.0` as a fallback to
454/// avoid division by zero.
455fn extract_diag_inv_f32(matrix: &CsrMatrix<f32>) -> Vec<f32> {
456    let n = matrix.rows;
457    let mut d_inv = vec![1.0f32; n];
458    for i in 0..n {
459        let start = matrix.row_ptr[i];
460        let end = matrix.row_ptr[i + 1];
461        for idx in start..end {
462            if matrix.col_indices[idx] == i {
463                let diag = matrix.values[idx];
464                if diag.abs() > 1e-15 {
465                    d_inv[i] = 1.0 / diag;
466                } else {
467                    warn!(
468                        row = i,
469                        diag_value = %diag,
470                        "zero or near-zero diagonal entry; substituting 1.0 — matrix may be singular"
471                    );
472                }
473                break;
474            }
475        }
476    }
477    d_inv
478}
479
480/// Compute the L2 (Euclidean) norm of a slice of `f32` values.
481///
482/// Uses `f64` accumulation to reduce catastrophic cancellation on large
483/// vectors.
484#[inline]
485fn l2_norm_f32(v: &[f32]) -> f32 {
486    let sum: f64 = v.iter().map(|&x| (x as f64) * (x as f64)).sum();
487    sum.sqrt() as f32
488}
489
490/// Scale every element of `v` by `s` in-place.
491#[inline]
492fn scale_vec_f32(v: &mut [f32], s: f32) {
493    for x in v.iter_mut() {
494        *x *= s;
495    }
496}
497
498// ---------------------------------------------------------------------------
499// Tests
500// ---------------------------------------------------------------------------
501
502#[cfg(test)]
503mod tests {
504    use super::*;
505    use crate::types::CsrMatrix;
506
507    /// Helper: build a diagonally dominant tridiagonal matrix.
508    fn tridiag_f32(n: usize, diag_val: f32, off_val: f32) -> CsrMatrix<f32> {
509        let mut entries = Vec::new();
510        for i in 0..n {
511            entries.push((i, i, diag_val));
512            if i > 0 {
513                entries.push((i, i - 1, off_val));
514            }
515            if i + 1 < n {
516                entries.push((i, i + 1, off_val));
517            }
518        }
519        CsrMatrix::<f32>::from_coo(n, n, entries)
520    }
521
522    /// Helper: build a 3x3 system whose eigenvalues are in (0, 2) so that
523    /// the Neumann series converges (rho(I - A) < 1).
524    fn test_matrix_f64() -> CsrMatrix<f64> {
525        CsrMatrix::<f64>::from_coo(
526            3,
527            3,
528            vec![
529                (0, 0, 1.0),
530                (0, 1, -0.1),
531                (1, 0, -0.1),
532                (1, 1, 1.0),
533                (1, 2, -0.1),
534                (2, 1, -0.1),
535                (2, 2, 1.0),
536            ],
537        )
538    }
539
540    #[test]
541    fn test_new() {
542        let solver = NeumannSolver::new(1e-8, 100);
543        assert_eq!(solver.tolerance, 1e-8);
544        assert_eq!(solver.max_iterations, 100);
545    }
546
547    #[test]
548    fn test_spectral_radius_identity() {
549        let identity = CsrMatrix::<f32>::identity(4);
550        let rho = NeumannSolver::estimate_spectral_radius(&identity);
551        assert!(rho < 0.1, "expected rho ~ 0 for identity, got {rho}");
552    }
553
554    #[test]
555    fn test_spectral_radius_pure_diagonal() {
556        // For a pure diagonal matrix D, D^{-1}A = I, so M = I - I = 0.
557        // The spectral radius should be ~0.
558        let a = CsrMatrix::<f32>::from_coo(3, 3, vec![(0, 0, 0.5_f32), (1, 1, 0.5), (2, 2, 0.5)]);
559        let rho = NeumannSolver::estimate_spectral_radius(&a);
560        assert!(rho < 0.1, "expected rho ~ 0 for diagonal matrix, got {rho}");
561    }
562
563    #[test]
564    fn test_spectral_radius_empty() {
565        let empty = CsrMatrix::<f32> {
566            row_ptr: vec![0],
567            col_indices: vec![],
568            values: vec![],
569            rows: 0,
570            cols: 0,
571        };
572        assert_eq!(NeumannSolver::estimate_spectral_radius(&empty), 0.0);
573    }
574
575    #[test]
576    fn test_spectral_radius_non_diag_dominant() {
577        // Matrix where off-diagonal entries dominate:
578        // [1  2]
579        // [2  1]
580        // D^{-1}A = [[1, 2], [2, 1]], so M = I - D^{-1}A = [[0, -2], [-2, 0]].
581        // Eigenvalues of M are +2 and -2, so rho(M) = 2 > 1.
582        let a = CsrMatrix::<f32>::from_coo(
583            2,
584            2,
585            vec![(0, 0, 1.0_f32), (0, 1, 2.0), (1, 0, 2.0), (1, 1, 1.0)],
586        );
587        let rho = NeumannSolver::estimate_spectral_radius(&a);
588        assert!(
589            rho > 1.0,
590            "expected rho > 1 for non-diag-dominant matrix, got {rho}"
591        );
592    }
593
594    #[test]
595    fn test_solve_identity() {
596        let identity = CsrMatrix::<f32>::identity(3);
597        let rhs = vec![1.0_f32, 2.0, 3.0];
598        let solver = NeumannSolver::new(1e-6, 100);
599        let result = solver.solve(&identity, &rhs).unwrap();
600        for (i, (&e, &a)) in rhs.iter().zip(result.solution.iter()).enumerate() {
601            assert!((e - a).abs() < 1e-4, "index {i}: expected {e}, got {a}");
602        }
603        assert!(result.residual_norm < 1e-6);
604    }
605
606    #[test]
607    fn test_solve_diagonal() {
608        let a = CsrMatrix::<f32>::from_coo(3, 3, vec![(0, 0, 0.5_f32), (1, 1, 0.5), (2, 2, 0.5)]);
609        let rhs = vec![1.0_f32, 1.0, 1.0];
610        let solver = NeumannSolver::new(1e-6, 200);
611        let result = solver.solve(&a, &rhs).unwrap();
612        for (i, &val) in result.solution.iter().enumerate() {
613            assert!(
614                (val - 2.0).abs() < 0.01,
615                "index {i}: expected ~2.0, got {val}"
616            );
617        }
618    }
619
620    #[test]
621    fn test_solve_tridiagonal() {
622        // diag=1.0, off=-0.1: Jacobi iteration matrix has rho ~ 0.17.
623        // Use 1e-6 tolerance since f32 accumulation limits floor.
624        let a = tridiag_f32(5, 1.0, -0.1);
625        let rhs = vec![1.0_f32, 0.0, 1.0, 0.0, 1.0];
626        let solver = NeumannSolver::new(1e-6, 1000);
627        let result = solver.solve(&a, &rhs).unwrap();
628        assert!(result.residual_norm < 1e-4);
629        assert!(result.iterations > 0);
630        assert!(!result.convergence_history.is_empty());
631    }
632
633    #[test]
634    fn test_solve_empty_system() {
635        let a = CsrMatrix::<f32> {
636            row_ptr: vec![0],
637            col_indices: vec![],
638            values: vec![],
639            rows: 0,
640            cols: 0,
641        };
642        let result = NeumannSolver::new(1e-6, 10).solve(&a, &[]).unwrap();
643        assert_eq!(result.iterations, 0);
644        assert!(result.solution.is_empty());
645    }
646
647    #[test]
648    fn test_solve_dimension_mismatch() {
649        let a = CsrMatrix::<f32>::identity(3);
650        let rhs = vec![1.0_f32, 2.0];
651        let err = NeumannSolver::new(1e-6, 100).solve(&a, &rhs).unwrap_err();
652        let msg = err.to_string();
653        assert!(
654            msg.contains("dimension") || msg.contains("mismatch"),
655            "got: {msg}"
656        );
657    }
658
659    #[test]
660    fn test_solve_non_square() {
661        let a = CsrMatrix::<f32>::from_coo(2, 3, vec![(0, 0, 1.0_f32), (1, 1, 1.0)]);
662        let rhs = vec![1.0_f32, 1.0];
663        let err = NeumannSolver::new(1e-6, 100).solve(&a, &rhs).unwrap_err();
664        let msg = err.to_string();
665        assert!(
666            msg.contains("square") || msg.contains("dimension"),
667            "got: {msg}"
668        );
669    }
670
671    #[test]
672    fn test_solve_divergent_matrix() {
673        // Non-diag-dominant: off-diagonal entries larger than diagonal.
674        let a = CsrMatrix::<f32>::from_coo(
675            2,
676            2,
677            vec![(0, 0, 1.0_f32), (0, 1, 2.0), (1, 0, 2.0), (1, 1, 1.0)],
678        );
679        let rhs = vec![1.0_f32, 1.0];
680        let err = NeumannSolver::new(1e-6, 100).solve(&a, &rhs).unwrap_err();
681        assert!(err.to_string().contains("spectral radius"), "got: {}", err);
682    }
683
684    #[test]
685    fn test_convergence_history_monotone() {
686        let a = CsrMatrix::<f32>::identity(4);
687        let rhs = vec![1.0_f32; 4];
688        let result = NeumannSolver::new(1e-10, 50).solve(&a, &rhs).unwrap();
689        assert!(!result.convergence_history.is_empty());
690        for window in result.convergence_history.windows(2) {
691            assert!(
692                window[1].residual_norm <= window[0].residual_norm + 1e-12,
693                "residual not decreasing: {} -> {}",
694                window[0].residual_norm,
695                window[1].residual_norm,
696            );
697        }
698    }
699
700    #[test]
701    fn test_algorithm_tag() {
702        let a = CsrMatrix::<f32>::identity(2);
703        let rhs = vec![1.0_f32; 2];
704        let result = NeumannSolver::new(1e-6, 100).solve(&a, &rhs).unwrap();
705        assert_eq!(result.algorithm, Algorithm::Neumann);
706    }
707
708    #[test]
709    fn test_solver_engine_trait_f64() {
710        let solver = NeumannSolver::new(1e-6, 200);
711        let engine: &dyn SolverEngine = &solver;
712        let a = test_matrix_f64();
713        let rhs = vec![1.0_f64, 0.0, 1.0];
714        let budget = ComputeBudget::default();
715        let result = engine.solve(&a, &rhs, &budget).unwrap();
716        assert!(result.residual_norm < 1e-4);
717        assert_eq!(result.algorithm, Algorithm::Neumann);
718    }
719
720    #[test]
721    fn test_larger_system_accuracy() {
722        let n = 50;
723        // diag=1.0, off=-0.1: Jacobi-preconditioned Neumann converges.
724        // Use 1e-6 tolerance for f32 precision headroom.
725        let a = tridiag_f32(n, 1.0, -0.1);
726        let rhs: Vec<f32> = (0..n).map(|i| (i as f32 + 1.0) / n as f32).collect();
727        let result = NeumannSolver::new(1e-6, 2000).solve(&a, &rhs).unwrap();
728        assert!(
729            result.residual_norm < 1e-6,
730            "residual too large: {}",
731            result.residual_norm
732        );
733        let mut ax = vec![0.0f32; n];
734        a.spmv(&result.solution, &mut ax);
735        for i in 0..n {
736            assert!(
737                (ax[i] - rhs[i]).abs() < 1e-4,
738                "A*x[{i}]={} but b[{i}]={}",
739                ax[i],
740                rhs[i]
741            );
742        }
743    }
744
745    #[test]
746    fn test_scalar_system() {
747        let a = CsrMatrix::<f32>::from_coo(1, 1, vec![(0, 0, 0.5_f32)]);
748        let rhs = vec![4.0_f32];
749        let result = NeumannSolver::new(1e-8, 200).solve(&a, &rhs).unwrap();
750        assert!(
751            (result.solution[0] - 8.0).abs() < 0.01,
752            "expected ~8.0, got {}",
753            result.solution[0]
754        );
755    }
756
757    #[test]
758    fn test_estimate_complexity() {
759        let solver = NeumannSolver::new(1e-6, 1000);
760        let profile = SparsityProfile {
761            rows: 100,
762            cols: 100,
763            nnz: 500,
764            density: 0.05,
765            is_diag_dominant: true,
766            estimated_spectral_radius: 0.5,
767            estimated_condition: 3.0,
768            is_symmetric_structure: true,
769            avg_nnz_per_row: 5.0,
770            max_nnz_per_row: 8,
771        };
772        let estimate = solver.estimate_complexity(&profile, 100);
773        assert_eq!(estimate.algorithm, Algorithm::Neumann);
774        assert!(estimate.estimated_flops > 0);
775        assert!(estimate.estimated_iterations > 0);
776        assert!(estimate.estimated_memory_bytes > 0);
777        assert_eq!(estimate.complexity_class, ComplexityClass::SublinearNnz);
778    }
779}