numeris 0.5.5

Pure-Rust numerical algorithms library — high performance with SIMD support while also supporting no-std for embedded and WASM targets.
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
use crate::linalg::{split_two_col_slices, LinalgError};
use crate::matrix::vector::Vector;
use crate::traits::{LinalgScalar, MatrixMut, MatrixRef};
use crate::Matrix;
use num_traits::Zero;

/// QR decomposition in place using Householder reflections.
///
/// On return, `a` contains the packed QR factorization:
/// - Upper triangle (including diagonal): R
/// - Lower triangle (excluding diagonal): Householder vectors (scaled)
///
/// `tau` is filled with the Householder scalar factors (length min(M,N)).
///
/// Works on rectangular matrices (M >= N).
/// Returns `LinalgError::Singular` if a zero column is encountered.
///
/// For complex matrices, uses `H = I - tau * v * v^H` (conjugate transpose).
pub fn qr_in_place<T: LinalgScalar>(
    a: &mut impl MatrixMut<T>,
    tau: &mut [T],
) -> Result<(), LinalgError> {
    let m = a.nrows();
    let n = a.ncols();
    let k = m.min(n);
    assert!(m >= n, "QR decomposition requires M >= N");
    assert_eq!(tau.len(), k, "tau length must equal min(M, N)");

    for col in 0..k {
        // Compute the squared norm of the sub-column a[col:m, col]
        // using |v|^2 = v * conj(v) for complex support
        let sub_col = a.col_as_slice(col, col);
        let mut norm_sq = <T::Real as Zero>::zero();
        for &v in sub_col {
            norm_sq = norm_sq + (v * v.conj()).re();
        }

        if norm_sq < T::lepsilon() {
            return Err(LinalgError::Singular);
        }

        let norm = norm_sq.lsqrt();
        let a_col_col = *a.get(col, col);

        // sigma: for real, sign(a[col,col]) * ||x||
        //        for complex, a[col,col] / |a[col,col]| * ||x||
        // This ensures v0 = a + sigma = (a/|a|)(|a| + ||x||) avoids cancellation.
        let alpha = a_col_col.modulus();
        let sigma = if alpha < T::lepsilon() {
            T::from_real(norm)
        } else {
            T::from_real(norm) * (a_col_col / T::from_real(alpha))
        };

        // v[col] = a[col,col] + sigma; rest of v is a[col+1:m, col] (stored in-place)
        let v0 = a_col_col + sigma;
        *a.get_mut(col, col) = v0;

        // tau = v0 / sigma
        // For real: tau = v0 / sigma (as before)
        // For complex: tau = 2 / ||v||^2 * v0^2 ... but the formula tau = v0/sigma
        // generalizes correctly when sigma = conj(a)/|a| * norm.
        let tau_val = v0 / sigma;
        tau[col] = tau_val;

        // Scale the sub-diagonal entries by 1/v0 for storage
        {
            let sub_col = a.col_as_mut_slice(col, col + 1);
            for x in sub_col.iter_mut() {
                *x = *x / v0;
            }
        }

        // Apply H to trailing columns: A[col:m, col+1:n] -= tau * v * (v^H * A)
        // where v = [1, a[col+1,col], ..., a[m-1,col]] (stored values)
        // v^H * A uses conj(v_i) for complex.
        //
        // Column-major: v sub-column and A[:,j] sub-column are contiguous.
        // The AXPY step (no conjugation) uses SIMD dispatch on the slices.
        for j in (col + 1)..n {
            // Compute v^H * A[:, j] over rows col..m (needs conj for complex)
            let mut dot = *a.get(col, j); // conj(1) * A[col,j] = A[col,j]
            let (v_slice, a_j_slice) = split_two_col_slices(a, col, j, col + 1);
            for idx in 0..v_slice.len() {
                dot = dot + v_slice[idx].conj() * a_j_slice[idx];
            }
            dot = dot * tau_val;

            // A[:, j] -= dot * v (no conjugation — uses SIMD AXPY dispatch)
            *a.get_mut(col, j) = *a.get(col, j) - dot; // v[0] = 1
            let (v_slice, a_j_slice) = split_two_col_slices(a, col, j, col + 1);
            crate::simd::axpy_neg_dispatch(a_j_slice, dot, v_slice);
        }

        // Store -sigma (the R diagonal entry) in a[col, col]
        *a.get_mut(col, col) = T::zero() - sigma;
    }

    Ok(())
}

/// QR decomposition of a fixed-size matrix (M >= N).
///
/// Stores the packed Householder vectors, R, and tau scalars.
/// Use `q()`, `r()`, `solve()`, or `det()` to work with the decomposition.
///
/// # Example
///
/// ```
/// use numeris::{Matrix, Vector};
///
/// // Least-squares fit: y = c0 + c1*x to points (0,1), (1,2), (2,4)
/// let a = Matrix::new([
///     [1.0_f64, 0.0],
///     [1.0, 1.0],
///     [1.0, 2.0],
/// ]);
/// let b = Vector::from_array([1.0, 2.0, 4.0]);
/// let x = a.qr().unwrap().solve(&b);
/// assert!((x[0] - 5.0 / 6.0).abs() < 1e-10);
/// assert!((x[1] - 3.0 / 2.0).abs() < 1e-10);
/// ```
#[derive(Debug)]
pub struct QrDecomposition<T, const M: usize, const N: usize> {
    qr: Matrix<T, M, N>,
    tau: [T; N],
}

impl<T: LinalgScalar, const M: usize, const N: usize> QrDecomposition<T, M, N> {
    /// Decompose a matrix. Returns an error if a column is rank-deficient.
    pub fn new(a: &Matrix<T, M, N>) -> Result<Self, LinalgError> {
        assert!(M >= N, "QR decomposition requires M >= N");
        let mut qr = *a;
        let mut tau = [T::zero(); N];
        qr_in_place(&mut qr, &mut tau)?;
        Ok(Self { qr, tau })
    }

    /// Extract the upper-triangular R factor (N × N).
    ///
    /// ```
    /// use numeris::Matrix;
    /// let a = Matrix::new([[12.0_f64, -51.0, 4.0], [6.0, 167.0, -68.0], [-4.0, 24.0, -41.0]]);
    /// let r = a.qr().unwrap().r();
    /// // R is upper-triangular
    /// assert!((r[(1, 0)]).abs() < 1e-12);
    /// assert!((r[(2, 0)]).abs() < 1e-12);
    /// ```
    pub fn r(&self) -> Matrix<T, N, N> {
        let mut r = Matrix::<T, N, N>::zeros();
        for i in 0..N {
            for j in i..N {
                r[(i, j)] = self.qr[(i, j)];
            }
        }
        r
    }

    /// Compute the thin Q factor (M × N, orthonormal columns).
    ///
    /// Applies Householder reflections in reverse to the first N columns
    /// of the identity matrix.
    ///
    /// For complex matrices, Q is unitary (`Q^H * Q = I`).
    ///
    /// ```
    /// use numeris::Matrix;
    /// let a = Matrix::new([[12.0_f64, -51.0, 4.0], [6.0, 167.0, -68.0], [-4.0, 24.0, -41.0]]);
    /// let qr = a.qr().unwrap();
    /// let q = qr.q();
    /// let qtq = q.transpose() * q;
    /// // Q^T * Q ≈ I
    /// assert!((qtq[(0, 0)] - 1.0).abs() < 1e-10);
    /// assert!((qtq[(0, 1)]).abs() < 1e-10);
    /// ```
    pub fn q(&self) -> Matrix<T, M, N> {
        // Start with the M×N "thin identity": e_0..e_{N-1}
        let mut q = Matrix::<T, M, N>::zeros();
        for i in 0..N {
            q[(i, i)] = T::one();
        }

        // Apply reflections in reverse order
        for col in (0..N).rev() {
            let tau_val = self.tau[col];

            // Apply H_col to Q[col:M, col:N]
            // v = [1, qr[col+1,col], ..., qr[M-1,col]]
            // H = I - tau * v * v^H
            let v_slice = self.qr.col_as_slice(col, col + 1);
            for j in col..N {
                let mut dot = q[(col, j)]; // conj(1) * q[col,j]
                let q_j_slice = q.col_as_slice(j, col + 1);
                for idx in 0..v_slice.len() {
                    dot = dot + v_slice[idx].conj() * q_j_slice[idx];
                }
                dot = dot * tau_val;

                q[(col, j)] = q[(col, j)] - dot;
                let q_j_slice = q.col_as_mut_slice(j, col + 1);
                crate::simd::axpy_neg_dispatch(q_j_slice, dot, v_slice);
            }
        }

        q
    }

    /// Solve `AX = B` for X where B has multiple columns (least-squares).
    ///
    /// Each column of B is solved independently.
    ///
    /// ```
    /// use numeris::Matrix;
    ///
    /// let a = Matrix::new([[1.0_f64, 0.0], [0.0, 1.0], [1.0, 1.0]]);
    /// let b = Matrix::new([[1.0, 0.0], [0.0, 1.0], [1.0, 1.0]]);
    /// let qr = a.qr().unwrap();
    /// let x = qr.solve_matrix(&b);
    /// assert_eq!(x.nrows(), 2);
    /// assert_eq!(x.ncols(), 2);
    /// ```
    pub fn solve_matrix<const P: usize>(&self, b: &Matrix<T, M, P>) -> Matrix<T, N, P> {
        let mut result = Matrix::<T, N, P>::zeros();
        for col in 0..P {
            let b_col = Vector::from_array(core::array::from_fn(|i| b[(i, col)]));
            let x_col = self.solve(&b_col);
            for i in 0..N {
                result[(i, col)] = x_col[i];
            }
        }
        result
    }

    /// Solve the least-squares problem min ||Ax - b|| for x.
    ///
    /// Computes x = R^{-1} Q^H b via Householder application + back substitution.
    pub fn solve(&self, b: &Vector<T, M>) -> Vector<T, N> {
        // Apply Q^H to b by applying each Householder reflection in order.
        // Work with a copy of b as a flat array of length M.
        let mut qtb = [T::zero(); M];
        for i in 0..M {
            qtb[i] = b[i];
        }

        for col in 0..N {
            let tau_val = self.tau[col];
            // v = [1, qr[col+1,col], ..., qr[M-1,col]]
            let v_slice = self.qr.col_as_slice(col, col + 1);
            let mut dot = qtb[col]; // conj(1) * qtb[col]
            for idx in 0..v_slice.len() {
                dot = dot + v_slice[idx].conj() * qtb[col + 1 + idx];
            }
            dot = dot * tau_val;

            qtb[col] = qtb[col] - dot;
            crate::simd::axpy_neg_dispatch(&mut qtb[col + 1..], dot, v_slice);
        }

        // Back substitution with R (upper triangle of qr, first N rows)
        let mut x = [T::zero(); N];
        for i in (0..N).rev() {
            let mut sum = qtb[i];
            for j in (i + 1)..N {
                sum = sum - self.qr[(i, j)] * x[j];
            }
            x[i] = sum / self.qr[(i, i)];
        }

        Vector::from_array(x)
    }

    /// Determinant of the original matrix (square only, M == N).
    ///
    /// Panics at runtime if M != N.
    pub fn det(&self) -> T {
        assert_eq!(M, N, "determinant requires a square matrix");
        let mut d = T::one();
        for i in 0..N {
            d = d * self.qr[(i, i)];
        }
        d
    }
}

// ── Column-pivoted QR ──────────────────────────────────────────────

/// QR decomposition with column pivoting (rank-revealing) in place.
///
/// On return, `a` contains the packed QR factorization (same layout as
/// `qr_in_place`), `tau` holds Householder scalars, and `perm[j]` gives the
/// original column index that was pivoted into column `j`.
///
/// At each step k, the column with the largest 2-norm in the remaining
/// sub-matrix is swapped into position k before applying the Householder
/// reflection.  Column norms are maintained incrementally and recomputed
/// when the incremental estimate drops below a threshold to avoid
/// cancellation drift.
///
/// Unlike `qr_in_place`, this function does **not** return an error for
/// rank-deficient matrices — it simply produces small/zero diagonal entries
/// in R, and the numerical rank can be read from those entries after the
/// fact.
pub fn qr_col_pivot_in_place<T: LinalgScalar>(
    a: &mut impl MatrixMut<T>,
    tau: &mut [T],
    perm: &mut [usize],
) {
    let m = a.nrows();
    let n = a.ncols();
    let k = m.min(n);
    assert!(m >= n, "QR decomposition requires M >= N");
    assert_eq!(tau.len(), k, "tau length must equal min(M, N)");
    assert_eq!(perm.len(), n, "perm length must equal N");

    // Initialise permutation as identity and column norms.
    // col_norms[j] tracks the squared 2-norm of the remaining sub-column
    // below the current pivot row.
    let mut col_norms = [<T::Real as Zero>::zero(); 64];
    // For matrices wider than 64 we fall back to a stack vec.  This is
    // an embedded-friendly library, so we avoid alloc here.
    assert!(n <= 64, "qr_col_pivot_in_place: N must be <= 64 for fixed-size stack storage");
    for j in 0..n {
        perm[j] = j;
        let col = a.col_as_slice(j, 0);
        let mut s = <T::Real as Zero>::zero();
        for &v in col {
            s = s + (v * v.conj()).re();
        }
        col_norms[j] = s;
    }

    for col in 0..k {
        // ── Pivot: find column with largest remaining norm ───────
        let mut best_j = col;
        let mut best_norm = col_norms[col];
        for j in (col + 1)..n {
            if col_norms[j] > best_norm {
                best_norm = col_norms[j];
                best_j = j;
            }
        }

        // Swap columns col <-> best_j in a, perm, and col_norms
        if best_j != col {
            perm.swap(col, best_j);
            col_norms.swap(col, best_j);
            // Swap entire columns in a
            for i in 0..m {
                let tmp = *a.get(i, col);
                *a.get_mut(i, col) = *a.get(i, best_j);
                *a.get_mut(i, best_j) = tmp;
            }
        }

        // ── Householder reflection on column col ────────────────
        let sub_col = a.col_as_slice(col, col);
        let mut norm_sq = <T::Real as Zero>::zero();
        for &v in sub_col {
            norm_sq = norm_sq + (v * v.conj()).re();
        }

        if norm_sq < T::lepsilon() {
            // Remaining columns are numerically zero — leave tau = 0
            tau[col] = T::zero();
            continue;
        }

        let norm = norm_sq.lsqrt();
        let a_col_col = *a.get(col, col);

        let alpha = a_col_col.modulus();
        let sigma = if alpha < T::lepsilon() {
            T::from_real(norm)
        } else {
            T::from_real(norm) * (a_col_col / T::from_real(alpha))
        };

        let v0 = a_col_col + sigma;
        *a.get_mut(col, col) = v0;

        let tau_val = v0 / sigma;
        tau[col] = tau_val;

        // Scale sub-diagonal entries
        {
            let sub_col = a.col_as_mut_slice(col, col + 1);
            for x in sub_col.iter_mut() {
                *x = *x / v0;
            }
        }

        // Apply H to trailing columns + update norms
        for j in (col + 1)..n {
            // v^H * A[:, j]
            let mut dot = *a.get(col, j);
            let (v_slice, a_j_slice) = split_two_col_slices(a, col, j, col + 1);
            for idx in 0..v_slice.len() {
                dot = dot + v_slice[idx].conj() * a_j_slice[idx];
            }
            dot = dot * tau_val;

            *a.get_mut(col, j) = *a.get(col, j) - dot;
            let (v_slice, a_j_slice) = split_two_col_slices(a, col, j, col + 1);
            crate::simd::axpy_neg_dispatch(a_j_slice, dot, v_slice);

            // Incremental norm update: subtract eliminated component
            let _eliminated = (*a.get(col, j) * (*a.get(col, j)).conj()).re();
            // This is the new a[col, j] after reflection, which became R[col, j].
            // The remaining norm below row col is col_norms[j] minus the
            // squared entry that was just "absorbed" into R row col.
            // But we need the norm of the sub-column below row col+1.
            // Recompute from scratch (safe, avoids cancellation issues):
            let sub = a.col_as_slice(j, col + 1);
            let mut s = <T::Real as Zero>::zero();
            for &v in sub {
                s = s + (v * v.conj()).re();
            }
            col_norms[j] = s;
        }

        // Store -sigma as R diagonal entry
        *a.get_mut(col, col) = T::zero() - sigma;
    }
}

/// Rank-revealing QR decomposition with column pivoting.
///
/// Factorises `A * P = Q * R` where P is a permutation matrix, Q is
/// orthogonal (unitary), and R is upper-triangular with diagonal entries
/// in decreasing magnitude.
///
/// The numerical rank can be determined by inspecting the R diagonal
/// via [`rank()`](QrPivotDecomposition::rank).
///
/// # Example
///
/// ```
/// use numeris::Matrix;
///
/// // Rank-2 matrix (col 2 = col 0 + col 1)
/// let a = Matrix::new([
///     [1.0_f64, 0.0, 1.0],
///     [0.0, 1.0, 1.0],
///     [1.0, 1.0, 2.0],
/// ]);
/// let qrp = a.qr_col_pivot();
/// assert_eq!(qrp.rank(1e-10), 2);
///
/// // Q * R reconstructs A with columns permuted
/// let q = qrp.q();
/// let r = qrp.r();
/// let qr = q * r;
/// let perm = qrp.permutation();
/// for i in 0..3 {
///     for j in 0..3 {
///         assert!((qr[(i, j)] - a[(i, perm[j])]).abs() < 1e-10);
///     }
/// }
/// ```
#[derive(Debug)]
pub struct QrPivotDecomposition<T, const M: usize, const N: usize> {
    qr: Matrix<T, M, N>,
    tau: [T; N],
    perm: [usize; N],
}

impl<T: LinalgScalar, const M: usize, const N: usize> QrPivotDecomposition<T, M, N> {
    /// Decompose a matrix with column pivoting.
    pub fn new(a: &Matrix<T, M, N>) -> Self {
        assert!(M >= N, "QR decomposition requires M >= N");
        let mut qr = *a;
        let mut tau = [T::zero(); N];
        let mut perm = [0usize; N];
        qr_col_pivot_in_place(&mut qr, &mut tau, &mut perm);
        Self { qr, tau, perm }
    }

    /// Extract the upper-triangular R factor (N x N).
    pub fn r(&self) -> Matrix<T, N, N> {
        let mut r = Matrix::<T, N, N>::zeros();
        for i in 0..N {
            for j in i..N {
                r[(i, j)] = self.qr[(i, j)];
            }
        }
        r
    }

    /// Compute the thin Q factor (M x N, orthonormal columns).
    pub fn q(&self) -> Matrix<T, M, N> {
        let mut q = Matrix::<T, M, N>::zeros();
        for i in 0..N {
            q[(i, i)] = T::one();
        }

        for col in (0..N).rev() {
            let tau_val = self.tau[col];
            if tau_val == T::zero() {
                continue;
            }

            let v_slice = self.qr.col_as_slice(col, col + 1);
            for j in col..N {
                let mut dot = q[(col, j)];
                let q_j_slice = q.col_as_slice(j, col + 1);
                for idx in 0..v_slice.len() {
                    dot = dot + v_slice[idx].conj() * q_j_slice[idx];
                }
                dot = dot * tau_val;

                q[(col, j)] = q[(col, j)] - dot;
                let q_j_slice = q.col_as_mut_slice(j, col + 1);
                crate::simd::axpy_neg_dispatch(q_j_slice, dot, v_slice);
            }
        }

        q
    }

    /// Column permutation vector.
    ///
    /// `perm[j]` is the index of the original column that was pivoted into
    /// position `j`.  The factorisation satisfies `A[:, perm] = Q * R`.
    pub fn permutation(&self) -> &[usize; N] {
        &self.perm
    }

    /// Numerical rank: number of R diagonal entries with magnitude above `tol`.
    pub fn rank(&self, tol: T::Real) -> usize {
        (0..N)
            .take_while(|&i| self.qr[(i, i)].modulus() > tol)
            .count()
    }
}

/// Convenience methods on rectangular matrices.
impl<T: LinalgScalar, const M: usize, const N: usize> Matrix<T, M, N> {
    /// QR decomposition using Householder reflections.
    ///
    /// Returns an error if the matrix is rank-deficient.
    pub fn qr(&self) -> Result<QrDecomposition<T, M, N>, LinalgError> {
        QrDecomposition::new(self)
    }

    /// Rank-revealing QR decomposition with column pivoting.
    ///
    /// Unlike [`qr()`](Matrix::qr), this never fails — rank-deficient matrices
    /// simply produce small diagonal entries in R.
    ///
    /// ```
    /// use numeris::Matrix;
    ///
    /// let a = Matrix::new([
    ///     [1.0_f64, 2.0, 3.0],
    ///     [4.0, 5.0, 6.0],
    ///     [7.0, 8.0, 9.0],
    /// ]);
    /// let qrp = a.qr_col_pivot();
    /// assert_eq!(qrp.rank(1e-10), 2); // rank-deficient
    /// ```
    pub fn qr_col_pivot(&self) -> QrPivotDecomposition<T, M, N> {
        QrPivotDecomposition::new(self)
    }
}

/// Convenience methods for QR solve on square matrices.
impl<T: LinalgScalar, const N: usize> Matrix<T, N, N> {
    /// Solve Ax = b via QR decomposition.
    pub fn solve_qr(&self, b: &Vector<T, N>) -> Result<Vector<T, N>, LinalgError> {
        Ok(self.qr()?.solve(b))
    }
}

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

    const TOL: f64 = 1e-10;

    fn assert_near(a: f64, b: f64, tol: f64, msg: &str) {
        assert!((a - b).abs() < tol, "{}: {} vs {} (diff {})", msg, a, b, (a - b).abs());
    }

    #[test]
    fn qr_square_3x3() {
        let a = Matrix::new([
            [12.0_f64, -51.0, 4.0],
            [6.0, 167.0, -68.0],
            [-4.0, 24.0, -41.0],
        ]);
        let qr = a.qr().unwrap();
        let q = qr.q();
        let r = qr.r();

        // Verify Q*R == A
        let qr_prod = q * r;
        for i in 0..3 {
            for j in 0..3 {
                assert_near(qr_prod[(i, j)], a[(i, j)], TOL,
                    &format!("QR[({},{})]", i, j));
            }
        }

        // Verify Q^T * Q == I
        let qtq = q.transpose() * q;
        for i in 0..3 {
            for j in 0..3 {
                let expected = if i == j { 1.0 } else { 0.0 };
                assert_near(qtq[(i, j)], expected, TOL,
                    &format!("QtQ[({},{})]", i, j));
            }
        }
    }

    #[test]
    fn qr_rectangular_4x3() {
        let a = Matrix::new([
            [1.0_f64, -1.0, 4.0],
            [1.0, 4.0, -2.0],
            [1.0, 4.0, 2.0],
            [1.0, -1.0, 0.0],
        ]);
        let qr = a.qr().unwrap();
        let q = qr.q();
        let r = qr.r();

        // Verify Q*R == A (Q is 4×3, R is 3×3)
        let qr_prod: Matrix<f64, 4, 3> = q * r;
        for i in 0..4 {
            for j in 0..3 {
                assert_near(qr_prod[(i, j)], a[(i, j)], TOL,
                    &format!("QR[({},{})]", i, j));
            }
        }

        // Verify Q^T * Q == I_3 (thin Q)
        let qtq: Matrix<f64, 3, 3> = q.transpose() * q;
        for i in 0..3 {
            for j in 0..3 {
                let expected = if i == j { 1.0 } else { 0.0 };
                assert_near(qtq[(i, j)], expected, TOL,
                    &format!("QtQ[({},{})]", i, j));
            }
        }
    }

    #[test]
    fn qr_solve_square() {
        // Compare QR solve against LU solve
        let a = Matrix::new([
            [2.0_f64, 1.0, -1.0],
            [-3.0, -1.0, 2.0],
            [-2.0, 1.0, 2.0],
        ]);
        let b = Vector::from_array([8.0, -11.0, -3.0]);

        let x_qr = a.solve_qr(&b).unwrap();
        let x_lu = a.solve(&b).unwrap();

        for i in 0..3 {
            assert_near(x_qr[i], x_lu[i], TOL, &format!("x[{}]", i));
        }
    }

    #[test]
    fn qr_least_squares() {
        // Overdetermined system: 3 equations, 2 unknowns
        // Fit y = c0 + c1*x to points (0,1), (1,2), (2,4)
        // A = [[1, 0], [1, 1], [1, 2]], b = [1, 2, 4]
        let a = Matrix::new([
            [1.0_f64, 0.0],
            [1.0, 1.0],
            [1.0, 2.0],
        ]);
        let b = Vector::from_array([1.0, 2.0, 4.0]);

        let qr = a.qr().unwrap();
        let x = qr.solve(&b);

        // Least-squares solution: A^T A x = A^T b
        // A^T A = [[3, 3], [3, 5]], A^T b = [7, 10]
        // x = [5/6, 3/2]
        assert_near(x[0], 5.0 / 6.0, TOL, "c0");
        assert_near(x[1], 3.0 / 2.0, TOL, "c1");

        // Verify the residual r = b - Ax is orthogonal to column space of A
        // i.e. A^T * r ≈ 0
        let ax = a * x;
        let r = b - ax;
        let at = a.transpose();
        let atr = at * r;
        for i in 0..2 {
            assert_near(atr[i], 0.0, TOL, &format!("A^T r[{}]", i));
        }
    }

    #[test]
    fn qr_det() {
        let a = Matrix::new([
            [6.0_f64, 1.0, 1.0],
            [4.0, -2.0, 5.0],
            [2.0, 8.0, 7.0],
        ]);
        let qr = a.qr().unwrap();
        let det_qr = qr.det();
        let det_lu = a.det();
        assert_near(det_qr.abs(), det_lu.abs(), TOL, "det magnitude");
    }

    #[test]
    fn qr_identity() {
        let id: Matrix<f64, 3, 3> = Matrix::eye();
        let qr = id.qr().unwrap();
        let q = qr.q();
        let r = qr.r();

        // Q should be identity (up to sign flips on columns)
        // R should be identity (up to sign flips on diagonal)
        // Q*R should be identity
        let prod = q * r;
        for i in 0..3 {
            for j in 0..3 {
                let expected = if i == j { 1.0 } else { 0.0 };
                assert_near(prod[(i, j)], expected, TOL,
                    &format!("QR[({},{})]", i, j));
            }
        }
    }

    #[test]
    fn qr_in_place_generic() {
        // Verify the free function works via the MatrixMut trait
        let mut a = Matrix::new([[2.0_f64, 1.0], [4.0, 3.0]]);
        let mut tau = [0.0; 2];
        let result = qr_in_place(&mut a, &mut tau);
        assert!(result.is_ok());
    }

    #[test]
    fn qr_rank_deficient() {
        // Matrix with a zero column
        let a = Matrix::new([
            [1.0_f64, 0.0],
            [0.0, 0.0],
        ]);
        assert_eq!(a.qr().unwrap_err(), LinalgError::Singular);
    }

    // ── Column-pivoted QR tests ────────────────────────────────────

    #[test]
    fn qr_pivot_full_rank_3x3() {
        let a = Matrix::new([
            [12.0_f64, -51.0, 4.0],
            [6.0, 167.0, -68.0],
            [-4.0, 24.0, -41.0],
        ]);
        let qrp = a.qr_col_pivot();
        let q = qrp.q();
        let r = qrp.r();
        let perm = qrp.permutation();

        // Q * R == A[:, perm]
        let qr_prod = q * r;
        for i in 0..3 {
            for j in 0..3 {
                assert_near(qr_prod[(i, j)], a[(i, perm[j])], TOL,
                    &format!("QR[({},{})]", i, j));
            }
        }

        // Q^T * Q == I
        let qtq = q.transpose() * q;
        for i in 0..3 {
            for j in 0..3 {
                let expected = if i == j { 1.0 } else { 0.0 };
                assert_near(qtq[(i, j)], expected, TOL,
                    &format!("QtQ[({},{})]", i, j));
            }
        }

        // Full rank
        assert_eq!(qrp.rank(1e-10), 3);
    }

    #[test]
    fn qr_pivot_rank_deficient() {
        // Rank-2 matrix: row 2 = row 0 + row 1
        let a = Matrix::new([
            [1.0_f64, 2.0, 3.0],
            [4.0, 5.0, 6.0],
            [5.0, 7.0, 9.0],
        ]);
        let qrp = a.qr_col_pivot();
        assert_eq!(qrp.rank(1e-10), 2);

        // Q * R still reconstructs A[:, perm]
        let q = qrp.q();
        let r = qrp.r();
        let perm = qrp.permutation();
        let qr_prod = q * r;
        for i in 0..3 {
            for j in 0..3 {
                assert_near(qr_prod[(i, j)], a[(i, perm[j])], TOL,
                    &format!("QR[({},{})]", i, j));
            }
        }
    }

    #[test]
    fn qr_pivot_rank_1() {
        // Rank-1 matrix: all rows are multiples of [1, 2, 3]
        let a = Matrix::new([
            [1.0_f64, 2.0, 3.0],
            [2.0, 4.0, 6.0],
            [3.0, 6.0, 9.0],
        ]);
        let qrp = a.qr_col_pivot();
        assert_eq!(qrp.rank(1e-10), 1);

        let q = qrp.q();
        let r = qrp.r();
        let perm = qrp.permutation();
        let qr_prod = q * r;
        for i in 0..3 {
            for j in 0..3 {
                assert_near(qr_prod[(i, j)], a[(i, perm[j])], TOL,
                    &format!("QR[({},{})]", i, j));
            }
        }
    }

    #[test]
    fn qr_pivot_zero_matrix() {
        let a = Matrix::<f64, 3, 3>::zeros();
        let qrp = a.qr_col_pivot();
        assert_eq!(qrp.rank(1e-10), 0);
    }

    #[test]
    fn qr_pivot_identity() {
        let id: Matrix<f64, 3, 3> = Matrix::eye();
        let qrp = id.qr_col_pivot();
        assert_eq!(qrp.rank(1e-10), 3);

        let q = qrp.q();
        let r = qrp.r();
        let perm = qrp.permutation();
        let qr_prod = q * r;
        for i in 0..3 {
            for j in 0..3 {
                assert_near(qr_prod[(i, j)], id[(i, perm[j])], TOL,
                    &format!("QR[({},{})]", i, j));
            }
        }
    }

    #[test]
    fn qr_pivot_rectangular_4x3() {
        let a = Matrix::new([
            [1.0_f64, -1.0, 4.0],
            [1.0, 4.0, -2.0],
            [1.0, 4.0, 2.0],
            [1.0, -1.0, 0.0],
        ]);
        let qrp = a.qr_col_pivot();
        let q = qrp.q();
        let r = qrp.r();
        let perm = qrp.permutation();

        // Q * R == A[:, perm]
        let qr_prod: Matrix<f64, 4, 3> = q * r;
        for i in 0..4 {
            for j in 0..3 {
                assert_near(qr_prod[(i, j)], a[(i, perm[j])], TOL,
                    &format!("QR[({},{})]", i, j));
            }
        }

        // Q^T * Q == I_3
        let qtq: Matrix<f64, 3, 3> = q.transpose() * q;
        for i in 0..3 {
            for j in 0..3 {
                let expected = if i == j { 1.0 } else { 0.0 };
                assert_near(qtq[(i, j)], expected, TOL,
                    &format!("QtQ[({},{})]", i, j));
            }
        }

        assert_eq!(qrp.rank(1e-10), 3);
    }

    #[test]
    fn qr_pivot_r_diagonal_decreasing() {
        // R diagonal magnitudes should be non-increasing
        let a = Matrix::new([
            [1.0_f64, 100.0, 0.5],
            [2.0, 200.0, 1.0],
            [3.0, 300.0, 1.5],
            [4.0, 400.0, 2.0],
        ]);
        let qrp = a.qr_col_pivot();
        let r = qrp.r();
        for i in 0..(3 - 1) {
            assert!(
                r[(i, i)].abs() >= r[(i + 1, i + 1)].abs() - TOL,
                "|R[{},{}]| = {} should >= |R[{},{}]| = {}",
                i, i,
                r[(i, i)].abs(),
                i + 1, i + 1,
                r[(i + 1, i + 1)].abs()
            );
        }
    }

    #[test]
    fn qr_pivot_permutation_is_valid() {
        let a = Matrix::new([
            [0.1_f64, 10.0, 5.0],
            [0.2, 20.0, 10.0],
            [0.3, 30.0, 15.0],
        ]);
        let qrp = a.qr_col_pivot();
        let perm = qrp.permutation();

        // permutation should contain each index exactly once
        let mut seen = [false; 3];
        for &p in perm {
            assert!(p < 3, "permutation index out of range");
            assert!(!seen[p], "duplicate permutation index");
            seen[p] = true;
        }
    }

    #[test]
    fn qr_pivot_column_dependent() {
        // Column 2 = column 0 + column 1
        let a = Matrix::new([
            [1.0_f64, 0.0, 1.0],
            [0.0, 1.0, 1.0],
            [1.0, 1.0, 2.0],
        ]);
        let qrp = a.qr_col_pivot();
        assert_eq!(qrp.rank(1e-10), 2);

        let q = qrp.q();
        let r = qrp.r();
        let perm = qrp.permutation();
        let qr_prod = q * r;
        for i in 0..3 {
            for j in 0..3 {
                assert_near(qr_prod[(i, j)], a[(i, perm[j])], TOL,
                    &format!("QR[({},{})]", i, j));
            }
        }
    }

    #[test]
    fn qr_pivot_2x2() {
        let a = Matrix::new([
            [3.0_f64, 1.0],
            [4.0, 2.0],
        ]);
        let qrp = a.qr_col_pivot();
        assert_eq!(qrp.rank(1e-10), 2);

        let q = qrp.q();
        let r = qrp.r();
        let perm = qrp.permutation();
        let qr_prod = q * r;
        for i in 0..2 {
            for j in 0..2 {
                assert_near(qr_prod[(i, j)], a[(i, perm[j])], TOL,
                    &format!("QR[({},{})]", i, j));
            }
        }
    }
}