ruvector-solver 2.0.6

Sublinear-time solver for RuVector: O(log n) to O(√n) algorithms for sparse linear systems, PageRank, and spectral methods
Documentation
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
//! Conjugate Gradient solver for symmetric positive-definite systems.
//!
//! Solves `Ax = b` where `A` is a symmetric positive-definite (SPD) sparse
//! matrix in CSR format. The algorithm converges in at most `n` iterations
//! for an `n x n` system, but in practice converges in
//! `O(sqrt(kappa) * log(1/eps))` iterations where `kappa = cond(A)`.
//!
//! # Algorithm
//!
//! Implements the Hestenes-Stiefel variant of Conjugate Gradient:
//!
//! ```text
//! r = b - A*x
//! z = M^{-1} * r        (preconditioner; z = r when disabled)
//! p = z
//! rz = r . z
//!
//! for k in 0..max_iterations:
//!     Ap = A * p
//!     alpha = rz / (p . Ap)
//!     x  = x + alpha * p
//!     r  = r - alpha * Ap
//!     if ||r||_2 < tolerance * ||b||_2:
//!         converged; break
//!     z  = M^{-1} * r
//!     rz_new = r . z
//!     beta = rz_new / rz
//!     p  = z + beta * p
//!     rz = rz_new
//! ```
//!
//! # Preconditioning
//!
//! When `use_preconditioner` is `true`, a diagonal (Jacobi) preconditioner is
//! applied: `M = diag(A)`, so that `z_i = r_i / A_{ii}`. This reduces the
//! effective condition number for diagonally-dominant systems without adding
//! significant per-iteration cost.
//!
//! # Numerical precision
//!
//! All dot products and norm computations use `f64` accumulation even though
//! the matrix may store `f32` values, preventing catastrophic cancellation in
//! the inner products that drive the CG recurrence.
//!
//! # Convergence
//!
//! Theoretical bound:
//! `||x_k - x*||_A <= 2 * ((sqrt(kappa) - 1)/(sqrt(kappa) + 1))^k * ||x_0 - x*||_A`
//! where `kappa` is the 2-condition number of `A`.

use std::time::Instant;

use tracing::{debug, trace, warn};

use crate::error::{SolverError, ValidationError};
use crate::traits::SolverEngine;
use crate::types::{
    Algorithm, ComplexityClass, ComplexityEstimate, ComputeBudget, ConvergenceInfo, CsrMatrix,
    SolverResult, SparsityProfile,
};

// ═══════════════════════════════════════════════════════════════════════════
// Helper functions -- f64-accumulated linear algebra primitives
// ═══════════════════════════════════════════════════════════════════════════

/// Compute the dot product of two `f32` slices using `f64` accumulation.
///
/// Uses a 4-wide accumulator to exploit instruction-level parallelism and
/// reduce the dependency chain length, preventing precision loss in the
/// inner products that drive the CG recurrence.
///
/// # Panics
///
/// Debug-asserts that `a.len() == b.len()`.
#[inline]
pub fn dot_product_f64(a: &[f32], b: &[f32]) -> f64 {
    assert_eq!(a.len(), b.len(), "dot_product_f64: length mismatch");

    let n = a.len();
    let chunks = n / 4;
    let remainder = n % 4;

    let mut acc0: f64 = 0.0;
    let mut acc1: f64 = 0.0;
    let mut acc2: f64 = 0.0;
    let mut acc3: f64 = 0.0;

    for i in 0..chunks {
        let j = i * 4;
        acc0 += a[j] as f64 * b[j] as f64;
        acc1 += a[j + 1] as f64 * b[j + 1] as f64;
        acc2 += a[j + 2] as f64 * b[j + 2] as f64;
        acc3 += a[j + 3] as f64 * b[j + 3] as f64;
    }

    let base = chunks * 4;
    for i in 0..remainder {
        acc0 += a[base + i] as f64 * b[base + i] as f64;
    }

    (acc0 + acc1) + (acc2 + acc3)
}

/// Compute the dot product of two `f64` slices with 4-wide accumulation.
///
/// Used internally when the working vectors are already `f64`.
#[inline]
fn dot_f64(a: &[f64], b: &[f64]) -> f64 {
    assert_eq!(a.len(), b.len(), "dot_f64: length mismatch");

    let n = a.len();
    let chunks = n / 4;
    let remainder = n % 4;

    let mut acc0: f64 = 0.0;
    let mut acc1: f64 = 0.0;
    let mut acc2: f64 = 0.0;
    let mut acc3: f64 = 0.0;

    for i in 0..chunks {
        let j = i * 4;
        acc0 += a[j] * b[j];
        acc1 += a[j + 1] * b[j + 1];
        acc2 += a[j + 2] * b[j + 2];
        acc3 += a[j + 3] * b[j + 3];
    }

    let base = chunks * 4;
    for i in 0..remainder {
        acc0 += a[base + i] * b[base + i];
    }

    (acc0 + acc1) + (acc2 + acc3)
}

/// Compute `y[i] += alpha * x[i]` for all `i` (AXPY operation, `f32`).
///
/// # Panics
///
/// Debug-asserts that `x.len() == y.len()`.
#[inline]
pub fn axpy(alpha: f32, x: &[f32], y: &mut [f32]) {
    assert_eq!(x.len(), y.len(), "axpy: length mismatch");

    let n = x.len();
    let chunks = n / 4;
    let base = chunks * 4;

    for i in 0..chunks {
        let j = i * 4;
        y[j] += alpha * x[j];
        y[j + 1] += alpha * x[j + 1];
        y[j + 2] += alpha * x[j + 2];
        y[j + 3] += alpha * x[j + 3];
    }
    for i in base..n {
        y[i] += alpha * x[i];
    }
}

/// Compute `y[i] += alpha * x[i]` for all `i` (AXPY operation, `f64`).
#[inline]
fn axpy_f64(alpha: f64, x: &[f64], y: &mut [f64]) {
    assert_eq!(x.len(), y.len(), "axpy_f64: length mismatch");

    let n = x.len();
    let chunks = n / 4;
    let base = chunks * 4;

    for i in 0..chunks {
        let j = i * 4;
        y[j] += alpha * x[j];
        y[j + 1] += alpha * x[j + 1];
        y[j + 2] += alpha * x[j + 2];
        y[j + 3] += alpha * x[j + 3];
    }
    for i in base..n {
        y[i] += alpha * x[i];
    }
}

/// Compute the L2 norm of an `f32` slice using `f64` accumulation.
///
/// Returns `sqrt(sum(x_i^2))` computed entirely in `f64`.
#[inline]
pub fn norm2(x: &[f32]) -> f64 {
    dot_product_f64(x, x).sqrt()
}

/// Compute the L2 norm of an `f64` slice.
#[inline]
fn norm2_f64(x: &[f64]) -> f64 {
    dot_f64(x, x).sqrt()
}

// ═══════════════════════════════════════════════════════════════════════════
// ConjugateGradientSolver
// ═══════════════════════════════════════════════════════════════════════════

/// Conjugate Gradient solver for symmetric positive-definite sparse systems.
///
/// Stores the solver configuration (tolerance, iteration cap, preconditioning).
/// The solve itself is stateless and may be invoked concurrently on different
/// inputs from multiple threads.
#[derive(Debug, Clone)]
pub struct ConjugateGradientSolver {
    /// Relative residual convergence tolerance.
    ///
    /// The solver declares convergence when `||r||_2 < tolerance * ||b||_2`.
    tolerance: f64,

    /// Maximum number of CG iterations before declaring non-convergence.
    max_iterations: usize,

    /// Whether to apply diagonal (Jacobi) preconditioning.
    ///
    /// When `true`, the preconditioner `M = diag(A)` is used, computing
    /// `z_i = r_i / A_{ii}` each iteration. Beneficial for diagonally-dominant
    /// systems at the cost of O(n) extra work per iteration.
    use_preconditioner: bool,
}

impl ConjugateGradientSolver {
    /// Create a new CG solver.
    ///
    /// # Arguments
    ///
    /// * `tolerance` -- Relative residual threshold for convergence. Must be
    ///   positive and finite.
    /// * `max_iterations` -- Upper bound on CG iterations. Must be >= 1.
    /// * `use_preconditioner` -- Enable diagonal (Jacobi) preconditioning.
    pub fn new(tolerance: f64, max_iterations: usize, use_preconditioner: bool) -> Self {
        Self {
            tolerance,
            max_iterations,
            use_preconditioner,
        }
    }

    /// Return the configured tolerance.
    #[inline]
    pub fn tolerance(&self) -> f64 {
        self.tolerance
    }

    /// Return the configured maximum iterations.
    #[inline]
    pub fn max_iterations(&self) -> usize {
        self.max_iterations
    }

    /// Return whether preconditioning is enabled.
    #[inline]
    pub fn use_preconditioner(&self) -> bool {
        self.use_preconditioner
    }

    // -------------------------------------------------------------------
    // Input validation
    // -------------------------------------------------------------------

    /// Validate inputs before entering the CG loop.
    fn validate(&self, matrix: &CsrMatrix<f64>, rhs: &[f64]) -> Result<(), SolverError> {
        if matrix.rows != matrix.cols {
            return Err(SolverError::InvalidInput(
                ValidationError::DimensionMismatch(format!(
                    "CG requires a square matrix but got {}x{}",
                    matrix.rows, matrix.cols,
                )),
            ));
        }

        if rhs.len() != matrix.rows {
            return Err(SolverError::InvalidInput(
                ValidationError::DimensionMismatch(format!(
                    "rhs length {} does not match matrix rows {}",
                    rhs.len(),
                    matrix.rows,
                )),
            ));
        }

        if matrix.row_ptr.len() != matrix.rows + 1 {
            return Err(SolverError::InvalidInput(
                ValidationError::DimensionMismatch(format!(
                    "row_ptr length {} does not equal rows + 1 = {}",
                    matrix.row_ptr.len(),
                    matrix.rows + 1,
                )),
            ));
        }

        if !self.tolerance.is_finite() || self.tolerance <= 0.0 {
            return Err(SolverError::InvalidInput(
                ValidationError::ParameterOutOfRange {
                    name: "tolerance".into(),
                    value: self.tolerance.to_string(),
                    expected: "positive finite value".into(),
                },
            ));
        }

        if self.max_iterations == 0 {
            return Err(SolverError::InvalidInput(
                ValidationError::ParameterOutOfRange {
                    name: "max_iterations".into(),
                    value: "0".into(),
                    expected: ">= 1".into(),
                },
            ));
        }

        Ok(())
    }

    // -------------------------------------------------------------------
    // Jacobi preconditioner
    // -------------------------------------------------------------------

    /// Build the Jacobi (diagonal) preconditioner from `A`.
    ///
    /// Returns `inv_diag[i] = 1.0 / A_{ii}`. Zero or near-zero diagonal
    /// entries are replaced with `1.0` to prevent division by zero.
    fn build_jacobi_preconditioner(matrix: &CsrMatrix<f64>) -> Vec<f64> {
        let n = matrix.rows;
        let mut inv_diag = vec![1.0f64; n];

        for row in 0..n {
            let start = matrix.row_ptr[row];
            let end = matrix.row_ptr[row + 1];
            for idx in start..end {
                if matrix.col_indices[idx] == row {
                    let diag_val = matrix.values[idx];
                    if diag_val.abs() > f64::EPSILON {
                        inv_diag[row] = 1.0 / diag_val;
                    }
                    break;
                }
            }
        }

        inv_diag
    }

    /// Apply the diagonal preconditioner: `z[i] = inv_diag[i] * r[i]`.
    #[inline]
    fn apply_preconditioner(inv_diag: &[f64], r: &[f64], z: &mut [f64]) {
        assert_eq!(inv_diag.len(), r.len());
        assert_eq!(r.len(), z.len());

        let n = r.len();
        let chunks = n / 4;
        let base = chunks * 4;

        for i in 0..chunks {
            let j = i * 4;
            z[j] = inv_diag[j] * r[j];
            z[j + 1] = inv_diag[j + 1] * r[j + 1];
            z[j + 2] = inv_diag[j + 2] * r[j + 2];
            z[j + 3] = inv_diag[j + 3] * r[j + 3];
        }
        for i in base..n {
            z[i] = inv_diag[i] * r[i];
        }
    }

    // -------------------------------------------------------------------
    // Core CG algorithm
    // -------------------------------------------------------------------

    /// Core CG algorithm implementation.
    ///
    /// Works entirely in `f64` precision internally. The final solution is
    /// down-cast to `f32` in the returned [`SolverResult`] to match the
    /// crate's output contract.
    fn solve_inner(
        &self,
        matrix: &CsrMatrix<f64>,
        rhs: &[f64],
        budget: &ComputeBudget,
    ) -> Result<SolverResult, SolverError> {
        let start_time = Instant::now();
        let n = matrix.rows;

        // Effective limits: take the tighter of our config and the budget.
        let effective_max_iter = self.max_iterations.min(budget.max_iterations);
        let effective_tol = self.tolerance.min(budget.tolerance);

        // --- Trivial case: zero-dimensional system ---
        if n == 0 {
            return Ok(SolverResult {
                solution: vec![],
                iterations: 0,
                residual_norm: 0.0,
                wall_time: start_time.elapsed(),
                convergence_history: vec![],
                algorithm: Algorithm::CG,
            });
        }

        // --- Allocate working vectors (all f64) ---
        let mut x = vec![0.0f64; n]; // solution (initial guess = zero)
        let mut r = vec![0.0f64; n]; // residual
        let mut z = vec![0.0f64; n]; // preconditioned residual
        let mut p = vec![0.0f64; n]; // search direction
        let mut ap = vec![0.0f64; n]; // A * p scratch buffer

        // --- Build preconditioner (if enabled) ---
        let inv_diag = if self.use_preconditioner {
            Some(Self::build_jacobi_preconditioner(matrix))
        } else {
            None
        };

        // --- r = b - A*x. Since x = 0, r = b ---
        r.copy_from_slice(rhs);

        // --- Convergence threshold: ||r||_2 < tol * ||b||_2 (relative) ---
        let b_norm = norm2_f64(rhs);
        let abs_tolerance = effective_tol * b_norm;

        // Handle zero RHS: the solution is the zero vector.
        if b_norm < f64::EPSILON {
            debug!("CG: zero RHS detected, returning zero solution");
            return Ok(SolverResult {
                solution: vec![0.0f32; n],
                iterations: 0,
                residual_norm: 0.0,
                wall_time: start_time.elapsed(),
                convergence_history: vec![],
                algorithm: Algorithm::CG,
            });
        }

        let initial_residual_norm = norm2_f64(&r);

        // --- z = M^{-1} * r ---
        match &inv_diag {
            Some(diag) => Self::apply_preconditioner(diag, &r, &mut z),
            None => z.copy_from_slice(&r),
        }

        // --- p = z ---
        p.copy_from_slice(&z);

        // --- rz = r . z ---
        let mut rz = dot_f64(&r, &z);

        let mut convergence_history = Vec::with_capacity(effective_max_iter.min(256));
        let mut converged = false;

        debug!(
            "CG: n={}, nnz={}, tol={:.2e}, max_iter={}, precond={}",
            n,
            matrix.nnz(),
            effective_tol,
            effective_max_iter,
            self.use_preconditioner,
        );

        // ===============================================================
        // Main CG loop (Hestenes-Stiefel)
        // ===============================================================
        for k in 0..effective_max_iter {
            // --- Budget: wall-time check ---
            if start_time.elapsed() > budget.max_time {
                warn!("CG: wall-time budget exhausted at iteration {k}");
                return Err(SolverError::BudgetExhausted {
                    reason: format!(
                        "wall-time limit {:?} exceeded at iteration {k}",
                        budget.max_time,
                    ),
                    elapsed: start_time.elapsed(),
                });
            }

            // --- Ap = A * p  (sparse matrix-vector product) ---
            matrix.spmv(&p, &mut ap);

            // --- alpha = rz / (p . Ap) ---
            let p_dot_ap = dot_f64(&p, &ap);

            // Guard: if p.Ap <= 0 the matrix is not SPD or we hit numerical
            // breakdown.
            if p_dot_ap <= 0.0 {
                warn!("CG: non-positive p.Ap = {p_dot_ap:.4e} at iteration {k}");
                return Err(SolverError::NumericalInstability {
                    iteration: k,
                    detail: format!("p.Ap = {p_dot_ap:.6e} <= 0; matrix may not be SPD",),
                });
            }

            let alpha = rz / p_dot_ap;

            // --- x = x + alpha * p ---
            axpy_f64(alpha, &p, &mut x);

            // --- r = r - alpha * Ap ---
            axpy_f64(-alpha, &ap, &mut r);

            // --- Convergence check: ||r||_2 < tol * ||b||_2 ---
            let r_norm = norm2_f64(&r);

            convergence_history.push(ConvergenceInfo {
                iteration: k,
                residual_norm: r_norm,
            });

            trace!(
                "CG iter {k}: ||r|| = {r_norm:.6e}, rel = {:.6e}",
                r_norm / b_norm,
            );

            if r_norm < abs_tolerance {
                converged = true;
                debug!(
                    "CG converged at iteration {k}: ||r|| = {r_norm:.6e}, \
                     rel = {:.6e}",
                    r_norm / b_norm,
                );
                break;
            }

            // --- Divergence detection ---
            // If ||r|| has grown by 10x from the initial residual, the
            // system is likely indefinite or ill-conditioned beyond rescue.
            if r_norm > 10.0 * initial_residual_norm {
                warn!(
                    "CG: divergence at iteration {k}: ||r|| = {r_norm:.6e} \
                     > 10 * ||r_0|| = {:.6e}",
                    10.0 * initial_residual_norm,
                );
                return Err(SolverError::NumericalInstability {
                    iteration: k,
                    detail: format!(
                        "residual diverged: ||r|| = {r_norm:.6e} exceeds \
                         10x initial residual {initial_residual_norm:.6e}",
                    ),
                });
            }

            // --- z = M^{-1} * r ---
            match &inv_diag {
                Some(diag) => Self::apply_preconditioner(diag, &r, &mut z),
                None => z.copy_from_slice(&r),
            }

            // --- rz_new = r . z ---
            let rz_new = dot_f64(&r, &z);

            // Guard: stagnation when rz is near-zero.
            if rz.abs() < f64::EPSILON * f64::EPSILON {
                warn!("CG: rz near zero at iteration {k}, stagnation");
                return Err(SolverError::NumericalInstability {
                    iteration: k,
                    detail: format!("rz = {rz:.6e} is near zero; solver stagnated",),
                });
            }

            // --- beta = rz_new / rz ---
            let beta = rz_new / rz;

            // --- p = z + beta * p ---
            for i in 0..n {
                p[i] = z[i] + beta * p[i];
            }

            // --- rz = rz_new ---
            rz = rz_new;
        }

        let wall_time = start_time.elapsed();
        let final_residual = norm2_f64(&r);

        if !converged {
            debug!(
                "CG: non-convergence after {} iterations, ||r|| = {final_residual:.6e}",
                effective_max_iter,
            );
            return Err(SolverError::NonConvergence {
                iterations: effective_max_iter,
                residual: final_residual,
                tolerance: abs_tolerance,
            });
        }

        // Down-cast the f64 solution to f32 for SolverResult.
        let solution_f32: Vec<f32> = x.iter().map(|&v| v as f32).collect();

        Ok(SolverResult {
            solution: solution_f32,
            iterations: convergence_history.len(),
            residual_norm: final_residual,
            wall_time,
            convergence_history,
            algorithm: Algorithm::CG,
        })
    }
}

// ═══════════════════════════════════════════════════════════════════════════
// SolverEngine trait implementation
// ═══════════════════════════════════════════════════════════════════════════

impl SolverEngine for ConjugateGradientSolver {
    /// Solve `Ax = b` using the Conjugate Gradient method.
    ///
    /// # Errors
    ///
    /// * [`SolverError::InvalidInput`] -- dimension mismatch or invalid params.
    /// * [`SolverError::NumericalInstability`] -- divergence or non-SPD matrix.
    /// * [`SolverError::NonConvergence`] -- iteration limit exceeded.
    /// * [`SolverError::BudgetExhausted`] -- wall-time limit exceeded.
    fn solve(
        &self,
        matrix: &CsrMatrix<f64>,
        rhs: &[f64],
        budget: &ComputeBudget,
    ) -> Result<SolverResult, SolverError> {
        self.validate(matrix, rhs)?;
        self.solve_inner(matrix, rhs, budget)
    }

    /// Estimate CG complexity from the sparsity profile.
    ///
    /// CG converges in `O(sqrt(kappa))` iterations, each costing `O(nnz)` for
    /// the SpMV plus `O(n)` for the vector updates.
    fn estimate_complexity(&self, profile: &SparsityProfile, n: usize) -> ComplexityEstimate {
        // Estimated iterations from condition number, clamped to max_iterations.
        let est_iters = (profile.estimated_condition.sqrt() as usize)
            .max(1)
            .min(self.max_iterations);

        // FLOPs per iteration: 2*nnz (SpMV) + 6*n (dot products, axpy ops).
        let flops_per_iter = 2 * profile.nnz as u64 + 6 * n as u64;
        let estimated_flops = est_iters as u64 * flops_per_iter;

        // Memory: 5 vectors of length n (x, r, z, p, Ap) plus preconditioner.
        let vec_bytes = n * std::mem::size_of::<f64>();
        let precond_bytes = if self.use_preconditioner {
            vec_bytes
        } else {
            0
        };
        let estimated_memory_bytes = 5 * vec_bytes + precond_bytes;

        ComplexityEstimate {
            algorithm: Algorithm::CG,
            estimated_flops,
            estimated_iterations: est_iters,
            estimated_memory_bytes,
            complexity_class: ComplexityClass::SqrtCondition,
        }
    }

    /// Return the algorithm identifier.
    fn algorithm(&self) -> Algorithm {
        Algorithm::CG
    }
}

// ═══════════════════════════════════════════════════════════════════════════
// Tests
// ═══════════════════════════════════════════════════════════════════════════

#[cfg(test)]
mod tests {
    use super::*;
    use std::time::Duration;

    /// Build a symmetric tridiagonal SPD matrix (f64):
    ///   diag = 4.0, off-diag = -1.0
    fn tridiagonal_spd(n: usize) -> CsrMatrix<f64> {
        let mut entries = Vec::with_capacity(3 * n);
        for i in 0..n {
            if i > 0 {
                entries.push((i, i - 1, -1.0f64));
            }
            entries.push((i, i, 4.0f64));
            if i + 1 < n {
                entries.push((i, i + 1, -1.0f64));
            }
        }
        CsrMatrix::<f64>::from_coo(n, n, entries)
    }

    /// Build a diagonal matrix from the given values.
    fn diagonal_matrix(diag: &[f64]) -> CsrMatrix<f64> {
        let n = diag.len();
        let entries: Vec<_> = diag.iter().enumerate().map(|(i, &v)| (i, i, v)).collect();
        CsrMatrix::<f64>::from_coo(n, n, entries)
    }

    /// Build an identity matrix of dimension n.
    fn identity(n: usize) -> CsrMatrix<f64> {
        CsrMatrix::<f64>::identity(n)
    }

    fn default_budget() -> ComputeBudget {
        ComputeBudget {
            max_time: Duration::from_secs(30),
            max_iterations: 10_000,
            tolerance: 1e-10,
        }
    }

    // -----------------------------------------------------------------
    // dot_product_f64
    // -----------------------------------------------------------------

    #[test]
    fn dot_product_f64_basic() {
        let a = vec![1.0f32, 2.0, 3.0];
        let b = vec![4.0f32, 5.0, 6.0];
        let result = dot_product_f64(&a, &b);
        assert!((result - 32.0).abs() < 1e-10);
    }

    #[test]
    fn dot_product_f64_empty() {
        assert!((dot_product_f64(&[], &[]) - 0.0).abs() < 1e-10);
    }

    #[test]
    fn dot_product_f64_precision() {
        let n = 10_000;
        let a = vec![1.0f32; n];
        let b = vec![1.0f32; n];
        assert!((dot_product_f64(&a, &b) - n as f64).abs() < 1e-10);
    }

    #[test]
    fn dot_product_f64_odd_length() {
        let a = vec![1.0f32, 2.0, 3.0, 4.0, 5.0];
        let b = vec![5.0f32, 4.0, 3.0, 2.0, 1.0];
        // 5 + 8 + 9 + 8 + 5 = 35
        assert!((dot_product_f64(&a, &b) - 35.0).abs() < 1e-10);
    }

    // -----------------------------------------------------------------
    // axpy
    // -----------------------------------------------------------------

    #[test]
    fn axpy_basic() {
        let x = vec![1.0f32, 2.0, 3.0];
        let mut y = vec![10.0f32, 20.0, 30.0];
        axpy(2.0, &x, &mut y);
        assert_eq!(y, vec![12.0, 24.0, 36.0]);
    }

    #[test]
    fn axpy_negative_alpha() {
        let x = vec![1.0f32, 1.0, 1.0];
        let mut y = vec![5.0f32, 5.0, 5.0];
        axpy(-3.0, &x, &mut y);
        assert_eq!(y, vec![2.0, 2.0, 2.0]);
    }

    // -----------------------------------------------------------------
    // norm2
    // -----------------------------------------------------------------

    #[test]
    fn norm2_basic() {
        let x = vec![3.0f32, 4.0];
        assert!((norm2(&x) - 5.0).abs() < 1e-10);
    }

    #[test]
    fn norm2_zero() {
        assert!((norm2(&vec![0.0f32; 5]) - 0.0).abs() < 1e-10);
    }

    // -----------------------------------------------------------------
    // CG solver: convergence on well-conditioned systems
    // -----------------------------------------------------------------

    #[test]
    fn cg_identity_matrix() {
        let n = 5;
        let matrix = identity(n);
        let rhs = vec![1.0, 2.0, 3.0, 4.0, 5.0];
        let budget = default_budget();

        let solver = ConjugateGradientSolver::new(1e-10, 100, false);
        let result = solver.solve(&matrix, &rhs, &budget).unwrap();

        for i in 0..n {
            assert!(
                (result.solution[i] as f64 - rhs[i]).abs() < 1e-5,
                "x[{i}] = {} != {}",
                result.solution[i],
                rhs[i],
            );
        }
        // Identity should converge in at most 1 iteration.
        assert!(result.iterations <= 1);
    }

    #[test]
    fn cg_diagonal_matrix() {
        let diag = vec![2.0, 3.0, 5.0, 7.0];
        let matrix = diagonal_matrix(&diag);
        let rhs = vec![4.0, 9.0, 25.0, 49.0];
        let budget = default_budget();

        let solver = ConjugateGradientSolver::new(1e-10, 100, false);
        let result = solver.solve(&matrix, &rhs, &budget).unwrap();

        let expected = [2.0, 3.0, 5.0, 7.0];
        for i in 0..4 {
            assert!(
                (result.solution[i] as f64 - expected[i]).abs() < 1e-4,
                "x[{i}] = {} != {}",
                result.solution[i],
                expected[i],
            );
        }
    }

    #[test]
    fn cg_tridiagonal_small() {
        let n = 10;
        let matrix = tridiagonal_spd(n);
        let rhs = vec![1.0f64; n];
        let budget = default_budget();

        let solver = ConjugateGradientSolver::new(1e-8, 200, false);
        let result = solver.solve(&matrix, &rhs, &budget).unwrap();

        assert!(
            result.residual_norm < 1e-6,
            "residual = {}",
            result.residual_norm,
        );
        assert!(
            result.iterations <= n,
            "took {} iterations for n={}",
            result.iterations,
            n,
        );
    }

    #[test]
    fn cg_tridiagonal_large() {
        let n = 500;
        let matrix = tridiagonal_spd(n);
        let rhs: Vec<f64> = (0..n).map(|i| (i as f64 + 1.0) / n as f64).collect();
        let budget = default_budget();

        let solver = ConjugateGradientSolver::new(1e-8, 2000, false);
        let result = solver.solve(&matrix, &rhs, &budget).unwrap();

        assert!(
            result.residual_norm < 1e-5,
            "residual = {}",
            result.residual_norm,
        );
    }

    // -----------------------------------------------------------------
    // Preconditioning
    // -----------------------------------------------------------------

    #[test]
    fn cg_preconditioned_converges_faster() {
        let n = 100;
        let matrix = tridiagonal_spd(n);
        let rhs = vec![1.0f64; n];
        let budget = default_budget();

        let no_precond = ConjugateGradientSolver::new(1e-8, 500, false);
        let with_precond = ConjugateGradientSolver::new(1e-8, 500, true);

        let result_no = no_precond.solve(&matrix, &rhs, &budget).unwrap();
        let result_yes = with_precond.solve(&matrix, &rhs, &budget).unwrap();

        assert!(result_no.residual_norm < 1e-6);
        assert!(result_yes.residual_norm < 1e-6);

        assert!(
            result_yes.iterations <= result_no.iterations,
            "preconditioned ({}) should use <= iterations than \
             unpreconditioned ({})",
            result_yes.iterations,
            result_no.iterations,
        );
    }

    // -----------------------------------------------------------------
    // Zero RHS / empty system
    // -----------------------------------------------------------------

    #[test]
    fn cg_zero_rhs() {
        let matrix = tridiagonal_spd(5);
        let rhs = vec![0.0f64; 5];
        let budget = default_budget();

        let solver = ConjugateGradientSolver::new(1e-8, 100, false);
        let result = solver.solve(&matrix, &rhs, &budget).unwrap();

        assert_eq!(result.iterations, 0);
        for &v in &result.solution {
            assert!((v as f64).abs() < 1e-10);
        }
    }

    #[test]
    fn cg_empty_system() {
        let matrix = CsrMatrix {
            row_ptr: vec![0],
            col_indices: vec![],
            values: Vec::<f64>::new(),
            rows: 0,
            cols: 0,
        };
        let rhs: Vec<f64> = vec![];
        let budget = default_budget();

        let solver = ConjugateGradientSolver::new(1e-8, 100, false);
        let result = solver.solve(&matrix, &rhs, &budget).unwrap();

        assert_eq!(result.iterations, 0);
        assert!(result.solution.is_empty());
    }

    // -----------------------------------------------------------------
    // Error cases
    // -----------------------------------------------------------------

    #[test]
    fn cg_dimension_mismatch() {
        let matrix = tridiagonal_spd(3);
        let rhs = vec![1.0f64; 5];
        let budget = default_budget();

        let solver = ConjugateGradientSolver::new(1e-8, 100, false);
        let err = solver.solve(&matrix, &rhs, &budget).unwrap_err();
        assert!(matches!(err, SolverError::InvalidInput(_)));
    }

    #[test]
    fn cg_non_square_matrix() {
        let matrix = CsrMatrix {
            row_ptr: vec![0, 1, 2],
            col_indices: vec![0, 1],
            values: vec![1.0f64, 1.0],
            rows: 2,
            cols: 3,
        };
        let rhs = vec![1.0f64; 2];
        let budget = default_budget();

        let solver = ConjugateGradientSolver::new(1e-8, 100, false);
        let err = solver.solve(&matrix, &rhs, &budget).unwrap_err();
        assert!(matches!(err, SolverError::InvalidInput(_)));
    }

    #[test]
    fn cg_non_convergence() {
        let n = 50;
        let matrix = tridiagonal_spd(n);
        let rhs = vec![1.0f64; n];
        let budget = ComputeBudget {
            max_time: Duration::from_secs(30),
            max_iterations: 1,
            tolerance: 1e-15,
        };

        let solver = ConjugateGradientSolver::new(1e-15, 1, false);
        let err = solver.solve(&matrix, &rhs, &budget).unwrap_err();
        assert!(matches!(err, SolverError::NonConvergence { .. }));
    }

    #[test]
    fn cg_budget_iteration_limit() {
        let n = 50;
        let matrix = tridiagonal_spd(n);
        let rhs = vec![1.0f64; n];

        // Solver allows 1000, but budget allows only 2.
        let solver = ConjugateGradientSolver::new(1e-15, 1000, false);
        let budget = ComputeBudget {
            max_time: Duration::from_secs(60),
            max_iterations: 2,
            tolerance: 1e-15,
        };

        let err = solver.solve(&matrix, &rhs, &budget).unwrap_err();
        assert!(matches!(err, SolverError::NonConvergence { .. }));
    }

    // -----------------------------------------------------------------
    // Convergence history
    // -----------------------------------------------------------------

    #[test]
    fn cg_convergence_history_populated() {
        let n = 20;
        let matrix = tridiagonal_spd(n);
        let rhs = vec![1.0f64; n];
        let budget = default_budget();

        let solver = ConjugateGradientSolver::new(1e-10, 200, false);
        let result = solver.solve(&matrix, &rhs, &budget).unwrap();

        assert!(!result.convergence_history.is_empty());

        // Final history entry should match the reported residual.
        let last = result.convergence_history.last().unwrap();
        assert!((last.residual_norm - result.residual_norm).abs() < 1e-12);
    }

    #[test]
    fn cg_algorithm_field() {
        let matrix = identity(3);
        let rhs = vec![1.0f64; 3];
        let budget = default_budget();

        let solver = ConjugateGradientSolver::new(1e-8, 100, false);
        let result = solver.solve(&matrix, &rhs, &budget).unwrap();
        assert_eq!(result.algorithm, Algorithm::CG);
    }

    // -----------------------------------------------------------------
    // Verify Ax = b for computed solution
    // -----------------------------------------------------------------

    #[test]
    fn cg_solution_satisfies_system() {
        let n = 20;
        let matrix = tridiagonal_spd(n);
        let rhs = vec![1.0f64; n];
        let budget = default_budget();

        let solver = ConjugateGradientSolver::new(1e-10, 200, true);
        let result = solver.solve(&matrix, &rhs, &budget).unwrap();

        // Up-cast solution back to f64 and compute Ax.
        let x_f64: Vec<f64> = result.solution.iter().map(|&v| v as f64).collect();
        let mut ax = vec![0.0f64; n];
        matrix.spmv(&x_f64, &mut ax);

        let mut max_err: f64 = 0.0;
        for i in 0..n {
            let err = (ax[i] - rhs[i]).abs();
            if err > max_err {
                max_err = err;
            }
        }

        assert!(
            max_err < 1e-4,
            "max |Ax - b| = {max_err:.6e}, expected < 1e-4",
        );
    }

    // -----------------------------------------------------------------
    // Estimate complexity
    // -----------------------------------------------------------------

    #[test]
    fn estimate_complexity_returns_cg() {
        let solver = ConjugateGradientSolver::new(1e-8, 500, true);
        let profile = SparsityProfile {
            rows: 100,
            cols: 100,
            nnz: 298,
            density: 0.0298,
            is_diag_dominant: true,
            estimated_spectral_radius: 0.5,
            estimated_condition: 100.0,
            is_symmetric_structure: true,
            avg_nnz_per_row: 2.98,
            max_nnz_per_row: 3,
        };

        let est = solver.estimate_complexity(&profile, 100);
        assert_eq!(est.algorithm, Algorithm::CG);
        assert_eq!(est.complexity_class, ComplexityClass::SqrtCondition);
        assert!(est.estimated_iterations > 0);
        assert!(est.estimated_flops > 0);
        assert!(est.estimated_memory_bytes > 0);
    }

    // -----------------------------------------------------------------
    // Accessors
    // -----------------------------------------------------------------

    #[test]
    fn accessors() {
        let solver = ConjugateGradientSolver::new(1e-6, 500, true);
        assert!((solver.tolerance() - 1e-6).abs() < 1e-15);
        assert_eq!(solver.max_iterations(), 500);
        assert!(solver.use_preconditioner());
    }
}