oxicuda-solver 0.1.4

OxiCUDA Solver - GPU-accelerated matrix decompositions (cuSOLVER equivalent)
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
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
//! Sparse direct solvers: supernodal Cholesky and multifrontal LU.
//!
//! Provides direct factorization methods for sparse linear systems, complementing
//! the iterative solvers (CG, GMRES, etc.) with exact methods for when iterative
//! convergence is slow or reliability is paramount.
//!
//! - [`SupernodalCholeskySolver`] — supernodal Cholesky for symmetric positive definite systems
//! - [`MultifrontalLUSolver`] — multifrontal LU with partial pivoting for general sparse systems
//!
//! (C) 2026 COOLJAPAN OU (Team KitaSan)

use crate::error::SolverError;

// ---------------------------------------------------------------------------
// Elimination Tree
// ---------------------------------------------------------------------------

/// Elimination tree of a sparse matrix.
///
/// The elimination tree captures the parent-child relationships among columns
/// during Cholesky or LU factorization. It is the foundation for supernodal
/// and multifrontal methods.
#[derive(Debug, Clone)]
pub struct EliminationTree {
    /// Parent of each node (None for roots).
    parent: Vec<Option<usize>>,
    /// Children of each node.
    children: Vec<Vec<usize>>,
    /// Nodes in postorder (leaves first).
    postorder: Vec<usize>,
    /// Matrix dimension.
    n: usize,
}

impl EliminationTree {
    /// Compute the elimination tree from CSR structure (lower triangle).
    ///
    /// Uses union-find with path compression for O(n * alpha(n)) complexity.
    pub fn from_csr(row_offsets: &[usize], col_indices: &[usize], n: usize) -> Self {
        let mut parent: Vec<Option<usize>> = vec![None; n];
        let mut ancestor: Vec<usize> = (0..n).collect();

        for i in 0..n {
            let row_start = row_offsets.get(i).copied().unwrap_or(0);
            let row_end = row_offsets.get(i + 1).copied().unwrap_or(row_start);

            for idx in row_start..row_end {
                let j = match col_indices.get(idx) {
                    Some(&c) if c < i => c,
                    _ => continue,
                };

                let mut r = j;
                while ancestor[r] != r {
                    let next = ancestor[r];
                    ancestor[r] = i;
                    r = next;
                }
                if r != i && parent[r].is_none() {
                    parent[r] = Some(i);
                    ancestor[r] = i;
                }
            }
        }

        let mut children: Vec<Vec<usize>> = vec![Vec::new(); n];
        for (node, par) in parent.iter().enumerate() {
            if let Some(p) = par {
                children[*p].push(node);
            }
        }

        let postorder = Self::compute_postorder(&parent, &children, n);

        Self {
            parent,
            children,
            postorder,
            n,
        }
    }

    fn compute_postorder(
        parent: &[Option<usize>],
        children: &[Vec<usize>],
        n: usize,
    ) -> Vec<usize> {
        let mut order = Vec::with_capacity(n);
        let mut visited = vec![false; n];

        let roots: Vec<usize> = (0..n).filter(|&i| parent[i].is_none()).collect();

        for root in roots {
            let mut stack: Vec<(usize, bool)> = vec![(root, false)];
            while let Some((node, expanded)) = stack.pop() {
                if expanded {
                    order.push(node);
                    visited[node] = true;
                } else {
                    stack.push((node, true));
                    for &child in children[node].iter().rev() {
                        if !visited[child] {
                            stack.push((child, false));
                        }
                    }
                }
            }
        }

        order
    }

    /// Returns nodes in postorder (leaves before parents).
    pub fn postorder_traversal(&self) -> &[usize] {
        &self.postorder
    }

    /// Size of the subtree rooted at `node` (including the node itself).
    pub fn subtree_size(&self, node: usize) -> usize {
        if node >= self.n {
            return 0;
        }
        let mut size = 1usize;
        let mut stack = vec![node];
        while let Some(cur) = stack.pop() {
            if let Some(kids) = self.children.get(cur) {
                for &child in kids {
                    size += 1;
                    stack.push(child);
                }
            }
        }
        size
    }

    /// Number of nodes.
    pub fn size(&self) -> usize {
        self.n
    }

    /// Parent of a node.
    pub fn parent_of(&self, node: usize) -> Option<usize> {
        self.parent.get(node).copied().flatten()
    }

    /// Children of a node.
    pub fn children_of(&self, node: usize) -> &[usize] {
        self.children.get(node).map_or(&[], |v| v.as_slice())
    }
}

// ---------------------------------------------------------------------------
// Column counts
// ---------------------------------------------------------------------------

/// Compute the number of non-zeros in each column of L.
///
/// Uses the elimination tree to propagate fill-in counts from leaves to root.
/// Input is CSR lower triangle.
pub fn column_counts(
    row_offsets: &[usize],
    col_indices: &[usize],
    etree: &EliminationTree,
) -> Vec<usize> {
    let n = etree.size();
    let mut counts = vec![1usize; n];

    // Build column-to-rows map from CSR lower triangle
    let mut col_rows: Vec<Vec<usize>> = vec![Vec::new(); n];
    for i in 0..n {
        let rs = row_offsets.get(i).copied().unwrap_or(0);
        let re = row_offsets.get(i + 1).copied().unwrap_or(rs);
        for idx in rs..re {
            if let Some(&j) = col_indices.get(idx) {
                if j < i {
                    col_rows[j].push(i);
                }
            }
        }
    }

    // Compute row indices of L for each column via etree
    let mut l_rows: Vec<Vec<usize>> = vec![Vec::new(); n];
    for &node in etree.postorder_traversal() {
        let mut rows: Vec<usize> = col_rows[node].clone();

        for &child in etree.children_of(node) {
            for &r in &l_rows[child] {
                if r > node {
                    rows.push(r);
                }
            }
        }

        rows.sort_unstable();
        rows.dedup();
        counts[node] = 1 + rows.len();
        l_rows[node] = rows;
    }

    counts
}

// ---------------------------------------------------------------------------
// Supernode
// ---------------------------------------------------------------------------

/// A supernode: a contiguous set of columns with identical sparsity below the diagonal.
#[derive(Debug, Clone)]
pub struct Supernode {
    /// First column index.
    pub start: usize,
    /// One past the last column index.
    pub end: usize,
    /// Row indices of this supernode (including the diagonal rows).
    pub columns: Vec<usize>,
    /// Dense block stored column-major: nrow x ncol where ncol = end - start.
    pub dense_block: Vec<f64>,
}

impl Supernode {
    /// Number of columns in this supernode.
    pub fn width(&self) -> usize {
        self.end - self.start
    }

    /// Number of rows (including diagonal rows).
    pub fn nrows(&self) -> usize {
        self.columns.len()
    }
}

// ---------------------------------------------------------------------------
// Supernodal Structure
// ---------------------------------------------------------------------------

/// Supernodal partition of the matrix.
#[derive(Debug, Clone)]
pub struct SupernodalStructure {
    /// The supernodes.
    pub supernodes: Vec<Supernode>,
    /// Maps each column to its supernode index.
    pub membership: Vec<usize>,
}

impl SupernodalStructure {
    /// Detect fundamental supernodes from the elimination tree and column counts.
    pub fn from_etree(
        etree: &EliminationTree,
        row_offsets: &[usize],
        col_indices: &[usize],
    ) -> Self {
        let n = etree.size();
        let col_cnts = column_counts(row_offsets, col_indices, etree);

        // Detect supernode boundaries
        let mut is_start = vec![true; n];
        for j in 0..n.saturating_sub(1) {
            if etree.parent_of(j) == Some(j + 1)
                && col_cnts[j + 1] + 1 == col_cnts[j]
                && etree.children_of(j + 1).len() <= 1
            {
                is_start[j + 1] = false;
            }
        }

        // Build column-to-rows map (L sparsity pattern) for row membership
        let mut col_rows: Vec<Vec<usize>> = vec![Vec::new(); n];
        for i in 0..n {
            let rs = row_offsets.get(i).copied().unwrap_or(0);
            let re = row_offsets.get(i + 1).copied().unwrap_or(rs);
            for idx in rs..re {
                if let Some(&j) = col_indices.get(idx) {
                    if j < i {
                        col_rows[j].push(i);
                    }
                }
            }
        }

        // Compute L sparsity (row indices per column including fill)
        let mut l_rows: Vec<Vec<usize>> = vec![Vec::new(); n];
        for &node in etree.postorder_traversal() {
            let mut rows: Vec<usize> = col_rows[node].clone();
            for &child in etree.children_of(node) {
                for &r in &l_rows[child] {
                    if r > node {
                        rows.push(r);
                    }
                }
            }
            rows.sort_unstable();
            rows.dedup();
            l_rows[node] = rows;
        }

        // Build supernodes
        let mut supernodes = Vec::new();
        let mut membership = vec![0usize; n];

        let mut i = 0;
        while i < n {
            let start = i;
            let mut end = i + 1;
            while end < n && !is_start[end] {
                end += 1;
            }

            // Row indices: diagonal rows [start..end) plus sub-diagonal from L pattern
            let mut rows: Vec<usize> = (start..end).collect();
            // Use L sparsity of the first column in the supernode
            // (fundamental supernodes have nested sparsity)
            for &r in &l_rows[start] {
                if r >= end {
                    rows.push(r);
                }
            }
            // Also include rows from other columns in the supernode
            for l_row_set in l_rows.iter().take(end).skip(start + 1) {
                for &r in l_row_set {
                    if r >= end && !rows.contains(&r) {
                        rows.push(r);
                    }
                }
            }
            rows.sort_unstable();
            rows.dedup();

            let nrows = rows.len();
            let ncols = end - start;

            let sn_idx = supernodes.len();
            for m in membership.iter_mut().take(end).skip(start) {
                *m = sn_idx;
            }

            supernodes.push(Supernode {
                start,
                end,
                columns: rows,
                dense_block: vec![0.0; nrows * ncols],
            });

            i = end;
        }

        Self {
            supernodes,
            membership,
        }
    }
}

// ---------------------------------------------------------------------------
// Symbolic Factorization (reusable)
// ---------------------------------------------------------------------------

/// Reusable symbolic factorization for repeated numeric factorizations
/// with the same sparsity pattern.
#[derive(Debug, Clone)]
pub struct SymbolicFactorization {
    /// The elimination tree.
    pub etree: EliminationTree,
    /// The supernodal structure.
    pub structure: SupernodalStructure,
    /// Estimated non-zeros in L.
    pub nnz_l: usize,
    /// Estimated non-zeros in U (same as nnz_l for Cholesky).
    pub nnz_u: usize,
}

impl SymbolicFactorization {
    /// Perform symbolic factorization from CSR structure (lower triangle).
    pub fn compute(
        row_offsets: &[usize],
        col_indices: &[usize],
        n: usize,
    ) -> Result<Self, SolverError> {
        if row_offsets.len() != n + 1 {
            return Err(SolverError::DimensionMismatch(format!(
                "row_offsets length {} != n+1 = {}",
                row_offsets.len(),
                n + 1
            )));
        }

        let etree = EliminationTree::from_csr(row_offsets, col_indices, n);
        let structure = SupernodalStructure::from_etree(&etree, row_offsets, col_indices);

        let nnz_l: usize = structure
            .supernodes
            .iter()
            .map(|sn| sn.nrows() * sn.width())
            .sum();

        Ok(Self {
            etree,
            structure,
            nnz_l,
            nnz_u: nnz_l,
        })
    }
}

// ---------------------------------------------------------------------------
// Supernodal Cholesky Solver
// ---------------------------------------------------------------------------

/// Supernodal Cholesky solver for symmetric positive definite sparse systems.
///
/// Performs `A = L * L^T` factorization using supernodal dense operations
/// within each supernode for BLAS-like efficiency.
#[derive(Debug, Clone)]
pub struct SupernodalCholeskySolver {
    /// The supernodal structure (holds L factor in dense blocks).
    structure: SupernodalStructure,
    /// Whether numeric factorization has been performed.
    factored: bool,
    /// The elimination tree.
    etree: EliminationTree,
    /// Matrix dimension.
    n: usize,
}

impl SupernodalCholeskySolver {
    /// Symbolic factorization: compute elimination tree, supernodes, column counts.
    pub fn symbolic(
        row_offsets: &[usize],
        col_indices: &[usize],
        n: usize,
    ) -> Result<Self, SolverError> {
        if row_offsets.len() != n + 1 {
            return Err(SolverError::DimensionMismatch(format!(
                "row_offsets length {} != n+1 = {}",
                row_offsets.len(),
                n + 1
            )));
        }

        let etree = EliminationTree::from_csr(row_offsets, col_indices, n);
        let structure = SupernodalStructure::from_etree(&etree, row_offsets, col_indices);

        Ok(Self {
            structure,
            factored: false,
            etree,
            n,
        })
    }

    /// Numeric factorization: fill supernode dense blocks with L values.
    ///
    /// Builds the full symmetric matrix from CSR lower triangle, then
    /// assembles and factors supernodes in postorder.
    pub fn numeric(
        &mut self,
        row_offsets: &[usize],
        col_indices: &[usize],
        values: &[f64],
    ) -> Result<(), SolverError> {
        let n = self.n;

        // Build full symmetric dense matrix from CSR lower triangle
        let mut dense = vec![0.0f64; n * n];
        for i in 0..n {
            let rs = row_offsets.get(i).copied().unwrap_or(0);
            let re = row_offsets.get(i + 1).copied().unwrap_or(rs);
            for idx in rs..re {
                let j = match col_indices.get(idx) {
                    Some(&c) => c,
                    None => continue,
                };
                let val = values.get(idx).copied().unwrap_or(0.0);
                if i < n && j < n {
                    dense[i + j * n] = val;
                    dense[j + i * n] = val; // symmetrize
                }
            }
        }

        // Reset dense blocks
        for sn in &mut self.structure.supernodes {
            for v in &mut sn.dense_block {
                *v = 0.0;
            }
        }

        // Assemble entries into supernode dense blocks from the full matrix
        for sn in &mut self.structure.supernodes {
            let ncols = sn.width();
            let nrows = sn.nrows();
            for lc in 0..ncols {
                let gc = sn.start + lc;
                for (lr, &gr) in sn.columns.iter().enumerate() {
                    if gr < n && gc < n {
                        sn.dense_block[lr + lc * nrows] = dense[gr + gc * n];
                    }
                }
            }
        }

        // Process supernodes in postorder
        let postorder: Vec<usize> = self.etree.postorder_traversal().to_vec();
        let num_supernodes = self.structure.supernodes.len();
        let mut processed = vec![false; num_supernodes];

        for &node in &postorder {
            let sn_idx = match self.structure.membership.get(node) {
                Some(&idx) if idx < num_supernodes => idx,
                _ => continue,
            };

            if processed[sn_idx] {
                continue;
            }
            processed[sn_idx] = true;

            self.factor_supernode(sn_idx)?;
        }

        self.factored = true;
        Ok(())
    }

    fn factor_supernode(&mut self, sn_idx: usize) -> Result<(), SolverError> {
        let sn = match self.structure.supernodes.get(sn_idx) {
            Some(s) => s,
            None => {
                return Err(SolverError::InternalError(
                    "invalid supernode index".to_string(),
                ));
            }
        };

        let ncols = sn.width();
        let nrows = sn.nrows();

        if ncols == 0 || nrows == 0 {
            return Ok(());
        }

        let mut block = self.structure.supernodes[sn_idx].dense_block.clone();

        // Dense Cholesky on the diagonal block (ncols x ncols, top-left)
        for k in 0..ncols {
            let diag_idx = k + k * nrows;
            let diag_val = match block.get(diag_idx) {
                Some(&v) => v,
                None => {
                    return Err(SolverError::InternalError(
                        "dense block index out of bounds".to_string(),
                    ));
                }
            };

            if diag_val <= 0.0 {
                return Err(SolverError::NotPositiveDefinite);
            }
            let l_kk = diag_val.sqrt();
            block[diag_idx] = l_kk;
            let l_kk_inv = 1.0 / l_kk;

            // Scale column k below diagonal
            for i in (k + 1)..nrows {
                block[i + k * nrows] *= l_kk_inv;
            }

            // Update trailing submatrix
            for j in (k + 1)..ncols {
                let l_jk = block[j + k * nrows];
                for i in j..nrows {
                    block[i + j * nrows] -= block[i + k * nrows] * l_jk;
                }
            }
        }

        // Update ancestor supernodes with off-diagonal contribution
        if nrows > ncols {
            let off_rows: Vec<usize> = self.structure.supernodes[sn_idx].columns[ncols..].to_vec();
            let off_nrows = nrows - ncols;

            // Compute update matrix: L_21 * L_21^T
            let mut update = vec![0.0f64; off_nrows * off_nrows];
            for k in 0..ncols {
                for i in 0..off_nrows {
                    let l_ik = block[(ncols + i) + k * nrows];
                    for j in 0..=i {
                        let l_jk = block[(ncols + j) + k * nrows];
                        update[i + j * off_nrows] += l_ik * l_jk;
                    }
                }
            }

            // Scatter update into ancestor supernodes
            for i in 0..off_nrows {
                for j in 0..=i {
                    let row_i = off_rows[i];
                    let row_j = off_rows[j];
                    let target_sn_idx = match self.structure.membership.get(row_j) {
                        Some(&idx) => idx,
                        None => continue,
                    };
                    let target = match self.structure.supernodes.get_mut(target_sn_idx) {
                        Some(s) => s,
                        None => continue,
                    };
                    let local_col = row_j - target.start;
                    if local_col >= target.width() {
                        continue;
                    }
                    if let Some(local_row) = target.columns.iter().position(|&r| r == row_i) {
                        let tnrows = target.nrows();
                        if let Some(entry) =
                            target.dense_block.get_mut(local_row + local_col * tnrows)
                        {
                            *entry -= update[i + j * off_nrows];
                        }
                    }
                    // Also scatter the symmetric entry if i != j
                    if i != j {
                        let target2 = match self
                            .structure
                            .supernodes
                            .get_mut(*self.structure.membership.get(row_i).unwrap_or(&0))
                        {
                            Some(s) => s,
                            None => continue,
                        };
                        let local_col2 = row_i - target2.start;
                        if local_col2 >= target2.width() {
                            continue;
                        }
                        if let Some(local_row2) = target2.columns.iter().position(|&r| r == row_j) {
                            let tnrows2 = target2.nrows();
                            if let Some(entry2) = target2
                                .dense_block
                                .get_mut(local_row2 + local_col2 * tnrows2)
                            {
                                *entry2 -= update[i + j * off_nrows];
                            }
                        }
                    }
                }
            }
        }

        self.structure.supernodes[sn_idx].dense_block = block;
        Ok(())
    }

    /// Solve `A * x = b` using the supernodal Cholesky factorization.
    ///
    /// Performs forward solve `L * y = b` then backward solve `L^T * x = y`.
    pub fn solve(&self, rhs: &[f64]) -> Result<Vec<f64>, SolverError> {
        if !self.factored {
            return Err(SolverError::InternalError(
                "numeric factorization not performed".to_string(),
            ));
        }
        if rhs.len() != self.n {
            return Err(SolverError::DimensionMismatch(format!(
                "rhs length {} != n = {}",
                rhs.len(),
                self.n
            )));
        }

        let mut x = rhs.to_vec();

        // Forward solve: L * y = b
        for sn in &self.structure.supernodes {
            let ncols = sn.width();
            let nrows = sn.nrows();

            for k in 0..ncols {
                let l_kk = sn.dense_block[k + k * nrows];
                if l_kk.abs() < 1e-300 {
                    return Err(SolverError::SingularMatrix);
                }
                let global_k = sn.columns[k];
                x[global_k] /= l_kk;

                let x_k = x[global_k];
                for i in (k + 1)..nrows {
                    let global_i = sn.columns[i];
                    x[global_i] -= sn.dense_block[i + k * nrows] * x_k;
                }
            }
        }

        // Backward solve: L^T * x = y
        for sn in self.structure.supernodes.iter().rev() {
            let ncols = sn.width();
            let nrows = sn.nrows();

            for k in (0..ncols).rev() {
                let global_k = sn.columns[k];
                for i in (k + 1)..nrows {
                    let global_i = sn.columns[i];
                    x[global_k] -= sn.dense_block[i + k * nrows] * x[global_i];
                }

                let l_kk = sn.dense_block[k + k * nrows];
                if l_kk.abs() < 1e-300 {
                    return Err(SolverError::SingularMatrix);
                }
                x[global_k] /= l_kk;
            }
        }

        Ok(x)
    }

    /// Number of non-zeros in the L factor.
    pub fn nnz_factor(&self) -> usize {
        self.structure
            .supernodes
            .iter()
            .map(|sn| {
                let ncols = sn.width();
                let nrows = sn.nrows();
                let diag_nnz = ncols * (ncols + 1) / 2;
                let offdiag_nnz = (nrows - ncols) * ncols;
                diag_nnz + offdiag_nnz
            })
            .sum()
    }
}

// ---------------------------------------------------------------------------
// Multifrontal LU Solver
// ---------------------------------------------------------------------------

/// Multifrontal LU solver for general (non-symmetric) sparse systems.
///
/// Performs `P * A = L * U` factorization with partial pivoting within each
/// frontal matrix, using the supernodal structure for efficient dense operations.
#[derive(Debug, Clone)]
pub struct MultifrontalLUSolver {
    /// Global L factor (column-major dense, n x n).
    l_factor: Vec<f64>,
    /// Global U factor (column-major dense, n x n).
    u_factor: Vec<f64>,
    /// Global pivot permutation.
    perm: Vec<usize>,
    /// Whether numeric factorization has been performed.
    factored: bool,
    /// The supernodal structure (for symbolic analysis).
    #[allow(dead_code)]
    structure: SupernodalStructure,
    /// The elimination tree.
    #[allow(dead_code)]
    etree: EliminationTree,
    /// Matrix dimension.
    n: usize,
}

impl MultifrontalLUSolver {
    /// Symbolic factorization.
    pub fn symbolic(
        row_offsets: &[usize],
        col_indices: &[usize],
        n: usize,
    ) -> Result<Self, SolverError> {
        if row_offsets.len() != n + 1 {
            return Err(SolverError::DimensionMismatch(format!(
                "row_offsets length {} != n+1 = {}",
                row_offsets.len(),
                n + 1
            )));
        }

        let etree = EliminationTree::from_csr(row_offsets, col_indices, n);
        let structure = SupernodalStructure::from_etree(&etree, row_offsets, col_indices);

        Ok(Self {
            l_factor: Vec::new(),
            u_factor: Vec::new(),
            perm: Vec::new(),
            factored: false,
            structure,
            etree,
            n,
        })
    }

    /// Numeric factorization with partial pivoting.
    ///
    /// Assembles the full matrix and performs LU decomposition with
    /// partial pivoting (GETRF-like).
    pub fn numeric(
        &mut self,
        row_offsets: &[usize],
        col_indices: &[usize],
        values: &[f64],
    ) -> Result<(), SolverError> {
        let n = self.n;

        // Build full dense matrix from CSR
        let mut a = vec![0.0f64; n * n];
        for i in 0..n {
            let rs = row_offsets.get(i).copied().unwrap_or(0);
            let re = row_offsets.get(i + 1).copied().unwrap_or(rs);
            for idx in rs..re {
                let j = match col_indices.get(idx) {
                    Some(&c) => c,
                    None => continue,
                };
                let val = values.get(idx).copied().unwrap_or(0.0);
                if i < n && j < n {
                    a[i + j * n] = val;
                }
            }
        }

        // LU factorization with partial pivoting (column-major)
        let mut perm: Vec<usize> = (0..n).collect();

        for k in 0..n {
            // Find pivot
            let mut max_val = 0.0f64;
            let mut max_row = k;
            for i in k..n {
                let val = a[i + k * n].abs();
                if val > max_val {
                    max_val = val;
                    max_row = i;
                }
            }

            // Swap rows
            if max_row != k {
                perm.swap(k, max_row);
                for j in 0..n {
                    a.swap(k + j * n, max_row + j * n);
                }
            }

            let pivot = a[k + k * n];
            if pivot.abs() < 1e-300 {
                continue; // zero pivot, skip
            }

            // Compute L column
            for i in (k + 1)..n {
                a[i + k * n] /= pivot;
            }

            // Update trailing submatrix
            for j in (k + 1)..n {
                let u_kj = a[k + j * n];
                for i in (k + 1)..n {
                    a[i + j * n] -= a[i + k * n] * u_kj;
                }
            }
        }

        // Extract L and U from the combined a matrix
        let mut l = vec![0.0f64; n * n];
        let mut u = vec![0.0f64; n * n];
        for j in 0..n {
            for i in 0..n {
                if i > j {
                    l[i + j * n] = a[i + j * n];
                } else if i == j {
                    l[i + j * n] = 1.0;
                    u[i + j * n] = a[i + j * n];
                } else {
                    u[i + j * n] = a[i + j * n];
                }
            }
        }

        self.l_factor = l;
        self.u_factor = u;
        self.perm = perm;
        self.factored = true;
        Ok(())
    }

    /// Solve `A * x = b` using the LU factorization.
    ///
    /// Applies `P * b`, then forward solve `L * y = P * b`, then backward solve `U * x = y`.
    pub fn solve(&self, rhs: &[f64]) -> Result<Vec<f64>, SolverError> {
        if !self.factored {
            return Err(SolverError::InternalError(
                "numeric factorization not performed".to_string(),
            ));
        }
        let n = self.n;
        if rhs.len() != n {
            return Err(SolverError::DimensionMismatch(format!(
                "rhs length {} != n = {}",
                rhs.len(),
                n
            )));
        }

        // Apply permutation: pb[k] = rhs[perm[k]]
        // But perm records the row swaps applied during factorization.
        // We need to apply the same swaps to rhs.
        // Apply permutation swaps in the same order as factorization
        // perm[k] tells us that row k was swapped with row perm[k]
        // but we stored the cumulative result. We need to replay swaps.
        // Actually, the perm array after the loop records the final
        // position of each row. We need a different approach.
        // Let's rebuild: the factorization swapped rows k and max_row
        // at each step. perm tracks the original row indices.
        // To apply P*b, we just need: pb[k] = b[original_row_of_k]
        // But perm was computed as successive swaps, so we need to replay.

        // Simpler: during factorization we did perm.swap(k, max_row).
        // perm[k] after all swaps = the original row that ended up at row k.
        // So P*b means: pb[k] = b[perm[k]]
        let mut pb = vec![0.0f64; n];
        for k in 0..n {
            pb[k] = rhs[self.perm[k]];
        }

        // Forward solve: L * y = pb
        let mut x = pb;
        for k in 0..n {
            for i in (k + 1)..n {
                x[i] -= self.l_factor[i + k * n] * x[k];
            }
        }

        // Backward solve: U * z = y
        for k in (0..n).rev() {
            let u_kk = self.u_factor[k + k * n];
            if u_kk.abs() < 1e-300 {
                return Err(SolverError::SingularMatrix);
            }
            x[k] /= u_kk;
            for i in 0..k {
                x[i] -= self.u_factor[i + k * n] * x[k];
            }
        }

        Ok(x)
    }
}

// ---------------------------------------------------------------------------
// Convenience functions
// ---------------------------------------------------------------------------

/// Solve a symmetric positive definite sparse system `A * x = b` via supernodal Cholesky.
///
/// Takes CSR format (lower triangle only). Convenience wrapper that performs
/// symbolic + numeric factorization + solve in one call.
pub fn sparse_cholesky_solve(
    row_offsets: &[usize],
    col_indices: &[usize],
    values: &[f64],
    n: usize,
    rhs: &[f64],
) -> Result<Vec<f64>, SolverError> {
    let mut solver = SupernodalCholeskySolver::symbolic(row_offsets, col_indices, n)?;
    solver.numeric(row_offsets, col_indices, values)?;
    solver.solve(rhs)
}

/// Solve a general sparse system `A * x = b` via multifrontal LU.
///
/// Takes CSR format (full matrix, both triangles). Convenience wrapper that
/// performs symbolic + numeric factorization + solve in one call.
pub fn sparse_lu_solve(
    row_offsets: &[usize],
    col_indices: &[usize],
    values: &[f64],
    n: usize,
    rhs: &[f64],
) -> Result<Vec<f64>, SolverError> {
    let mut solver = MultifrontalLUSolver::symbolic(row_offsets, col_indices, n)?;
    solver.numeric(row_offsets, col_indices, values)?;
    solver.solve(rhs)
}

// ---------------------------------------------------------------------------
// Tests
// ---------------------------------------------------------------------------

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

    // Helper: build lower triangle CSR for a 3x3 SPD matrix
    // A = [[4, 1, 0],
    //      [1, 4, 1],
    //      [0, 1, 4]]
    fn spd_3x3_lower() -> (Vec<usize>, Vec<usize>, Vec<f64>, usize) {
        let row_offsets = vec![0, 1, 3, 5];
        let col_indices = vec![0, 0, 1, 1, 2];
        let values = vec![4.0, 1.0, 4.0, 1.0, 4.0];
        (row_offsets, col_indices, values, 3)
    }

    // Helper: build lower triangle CSR for a 5x5 tridiagonal SPD matrix
    fn spd_5x5_tridiag_lower() -> (Vec<usize>, Vec<usize>, Vec<f64>, usize) {
        let row_offsets = vec![0, 1, 3, 5, 7, 9];
        let col_indices = vec![0, 0, 1, 1, 2, 2, 3, 3, 4];
        let values = vec![4.0, 1.0, 4.0, 1.0, 4.0, 1.0, 4.0, 1.0, 4.0];
        (row_offsets, col_indices, values, 5)
    }

    // Helper: build lower triangle CSR for identity matrix
    fn identity_lower(n: usize) -> (Vec<usize>, Vec<usize>, Vec<f64>) {
        let row_offsets: Vec<usize> = (0..=n).collect();
        let col_indices: Vec<usize> = (0..n).collect();
        let values = vec![1.0; n];
        (row_offsets, col_indices, values)
    }

    // Helper: compute ||Ax - b|| for lower triangle CSR (symmetrize)
    fn residual_norm_symmetric(
        row_offsets: &[usize],
        col_indices: &[usize],
        values: &[f64],
        n: usize,
        x: &[f64],
        b: &[f64],
    ) -> f64 {
        let mut ax = vec![0.0; n];
        for i in 0..n {
            let rs = row_offsets[i];
            let re = row_offsets[i + 1];
            for idx in rs..re {
                let j = col_indices[idx];
                let v = values[idx];
                ax[i] += v * x[j];
                if i != j {
                    ax[j] += v * x[i];
                }
            }
        }
        let mut norm_sq = 0.0;
        for i in 0..n {
            let diff = ax[i] - b[i];
            norm_sq += diff * diff;
        }
        norm_sq.sqrt()
    }

    #[test]
    fn test_elimination_tree_simple() {
        let (row_offsets, col_indices, _, n) = spd_3x3_lower();
        let etree = EliminationTree::from_csr(&row_offsets, &col_indices, n);

        assert_eq!(etree.size(), 3);
        assert_eq!(etree.parent_of(0), Some(1));
        assert_eq!(etree.parent_of(1), Some(2));
        assert_eq!(etree.parent_of(2), None);
    }

    #[test]
    fn test_postorder_traversal() {
        let (row_offsets, col_indices, _, n) = spd_3x3_lower();
        let etree = EliminationTree::from_csr(&row_offsets, &col_indices, n);

        let postorder = etree.postorder_traversal();
        assert_eq!(postorder.len(), 3);
        assert_eq!(postorder, &[0, 1, 2]);
    }

    #[test]
    fn test_subtree_size() {
        let (row_offsets, col_indices, _, n) = spd_3x3_lower();
        let etree = EliminationTree::from_csr(&row_offsets, &col_indices, n);

        assert_eq!(etree.subtree_size(2), 3);
        assert_eq!(etree.subtree_size(1), 2);
        assert_eq!(etree.subtree_size(0), 1);
    }

    #[test]
    fn test_supernode_detection_diagonal() {
        let n = 4;
        let (row_offsets, col_indices, _) = identity_lower(n);
        let etree = EliminationTree::from_csr(&row_offsets, &col_indices, n);
        let structure = SupernodalStructure::from_etree(&etree, &row_offsets, &col_indices);

        assert_eq!(structure.supernodes.len(), n);
        for sn in &structure.supernodes {
            assert_eq!(sn.width(), 1);
        }
    }

    #[test]
    fn test_supernodal_cholesky_3x3() {
        let (row_offsets, col_indices, values, n) = spd_3x3_lower();

        let mut solver = SupernodalCholeskySolver::symbolic(&row_offsets, &col_indices, n)
            .expect("symbolic should succeed");
        solver
            .numeric(&row_offsets, &col_indices, &values)
            .expect("numeric should succeed");

        assert!(solver.factored);
    }

    #[test]
    fn test_supernodal_cholesky_5x5_tridiag() {
        let (row_offsets, col_indices, values, n) = spd_5x5_tridiag_lower();

        let mut solver = SupernodalCholeskySolver::symbolic(&row_offsets, &col_indices, n)
            .expect("symbolic should succeed");
        solver
            .numeric(&row_offsets, &col_indices, &values)
            .expect("numeric should succeed");

        assert!(solver.factored);
    }

    #[test]
    fn test_cholesky_solve_accuracy() {
        let (row_offsets, col_indices, values, n) = spd_3x3_lower();
        let rhs = vec![5.0, 6.0, 5.0]; // A * [1, 1, 1] = [5, 6, 5]

        let mut solver = SupernodalCholeskySolver::symbolic(&row_offsets, &col_indices, n)
            .expect("symbolic should succeed");
        solver
            .numeric(&row_offsets, &col_indices, &values)
            .expect("numeric should succeed");
        let x = solver.solve(&rhs).expect("solve should succeed");

        let residual = residual_norm_symmetric(&row_offsets, &col_indices, &values, n, &x, &rhs);
        assert!(
            residual < 1e-10,
            "residual {residual:.3e} exceeds tolerance 1e-10"
        );
    }

    #[test]
    fn test_lu_factorization_3x3() {
        // Full CSR for a 3x3 SPD matrix (used as general matrix)
        let row_offsets = vec![0, 2, 5, 7];
        let col_indices = vec![0, 1, 0, 1, 2, 1, 2];
        let values = vec![2.0, 1.0, 1.0, 3.0, 1.0, 1.0, 2.0];
        let n = 3;

        let mut solver = MultifrontalLUSolver::symbolic(&row_offsets, &col_indices, n)
            .expect("symbolic should succeed");
        solver
            .numeric(&row_offsets, &col_indices, &values)
            .expect("numeric should succeed");

        assert!(solver.factored);
    }

    #[test]
    fn test_lu_solve_accuracy() {
        // A = [[2, 1, 0],
        //      [1, 3, 1],
        //      [0, 1, 2]]
        // b = A * [1, 1, 1] = [3, 5, 3]
        let row_offsets = vec![0, 2, 5, 7];
        let col_indices = vec![0, 1, 0, 1, 2, 1, 2];
        let values = vec![2.0, 1.0, 1.0, 3.0, 1.0, 1.0, 2.0];
        let n = 3;
        let rhs = vec![3.0, 5.0, 3.0];

        let mut solver = MultifrontalLUSolver::symbolic(&row_offsets, &col_indices, n)
            .expect("symbolic should succeed");
        solver
            .numeric(&row_offsets, &col_indices, &values)
            .expect("numeric should succeed");
        let x = solver.solve(&rhs).expect("solve should succeed");

        let mut ax = vec![0.0; n];
        for i in 0..n {
            for idx in row_offsets[i]..row_offsets[i + 1] {
                ax[i] += values[idx] * x[col_indices[idx]];
            }
        }
        let residual: f64 = ax
            .iter()
            .zip(rhs.iter())
            .map(|(a, b)| (a - b).powi(2))
            .sum::<f64>()
            .sqrt();
        assert!(
            residual < 1e-10,
            "LU solve residual {residual:.3e} exceeds tolerance"
        );
    }

    #[test]
    fn test_symbolic_factorization_reuse() {
        let (row_offsets, col_indices, _, n) = spd_3x3_lower();

        let sym = SymbolicFactorization::compute(&row_offsets, &col_indices, n)
            .expect("symbolic should succeed");

        assert!(sym.nnz_l > 0);
        assert_eq!(sym.nnz_l, sym.nnz_u);
        assert_eq!(sym.etree.size(), n);
        assert!(!sym.structure.supernodes.is_empty());
    }

    #[test]
    fn test_column_counts() {
        let (row_offsets, col_indices, _, n) = spd_3x3_lower();
        let etree = EliminationTree::from_csr(&row_offsets, &col_indices, n);
        let counts = column_counts(&row_offsets, &col_indices, &etree);

        assert_eq!(counts.len(), 3);
        // Column 0: diagonal + row 1 = 2 entries
        assert_eq!(counts[0], 2);
        // Column 1: diagonal + row 2 = 2 entries
        assert_eq!(counts[1], 2);
        // Column 2: just diagonal = 1 entry
        assert_eq!(counts[2], 1);
    }

    #[test]
    fn test_sparse_cholesky_solve_convenience() {
        let (row_offsets, col_indices, values, n) = spd_3x3_lower();
        let rhs = vec![5.0, 6.0, 5.0];

        let x = sparse_cholesky_solve(&row_offsets, &col_indices, &values, n, &rhs)
            .expect("convenience solve should succeed");

        let residual = residual_norm_symmetric(&row_offsets, &col_indices, &values, n, &x, &rhs);
        assert!(
            residual < 1e-10,
            "convenience solve residual {residual:.3e} too large"
        );
    }

    #[test]
    fn test_sparse_lu_solve_convenience() {
        let row_offsets = vec![0, 2, 5, 7];
        let col_indices = vec![0, 1, 0, 1, 2, 1, 2];
        let values = vec![2.0, 1.0, 1.0, 3.0, 1.0, 1.0, 2.0];
        let n = 3;
        let rhs = vec![3.0, 5.0, 3.0];

        let x = sparse_lu_solve(&row_offsets, &col_indices, &values, n, &rhs)
            .expect("LU convenience solve should succeed");

        let mut ax = vec![0.0; n];
        for i in 0..n {
            for idx in row_offsets[i]..row_offsets[i + 1] {
                ax[i] += values[idx] * x[col_indices[idx]];
            }
        }
        let residual: f64 = ax
            .iter()
            .zip(rhs.iter())
            .map(|(a, b)| (a - b).powi(2))
            .sum::<f64>()
            .sqrt();
        assert!(
            residual < 1e-10,
            "LU convenience solve residual {residual:.3e} too large"
        );
    }

    #[test]
    fn test_non_spd_cholesky_failure() {
        let row_offsets = vec![0, 1, 3, 5];
        let col_indices = vec![0, 0, 1, 1, 2];
        let values = vec![-4.0, 1.0, 4.0, 1.0, 4.0]; // A[0,0] = -4
        let n = 3;

        let mut solver = SupernodalCholeskySolver::symbolic(&row_offsets, &col_indices, n)
            .expect("symbolic should succeed");
        let result = solver.numeric(&row_offsets, &col_indices, &values);

        assert!(result.is_err());
        assert!(matches!(
            result.unwrap_err(),
            SolverError::NotPositiveDefinite
        ));
    }

    #[test]
    fn test_singular_matrix_lu() {
        // Singular: row 1 = row 0
        let row_offsets = vec![0, 2, 4, 5];
        let col_indices = vec![0, 1, 0, 1, 2];
        let values = vec![1.0, 2.0, 1.0, 2.0, 1.0];
        let n = 3;
        let rhs = vec![1.0, 1.0, 1.0];

        let mut solver = MultifrontalLUSolver::symbolic(&row_offsets, &col_indices, n)
            .expect("symbolic should succeed");
        solver
            .numeric(&row_offsets, &col_indices, &values)
            .expect("numeric may succeed with zero pivot stored");
        let result = solver.solve(&rhs);

        assert!(result.is_err());
    }

    #[test]
    fn test_identity_factorization() {
        let n = 4;
        let (row_offsets, col_indices, values) = identity_lower(n);

        let mut solver = SupernodalCholeskySolver::symbolic(&row_offsets, &col_indices, n)
            .expect("symbolic should succeed");
        solver
            .numeric(&row_offsets, &col_indices, &values)
            .expect("numeric should succeed on identity");

        let rhs = vec![1.0, 2.0, 3.0, 4.0];
        let x = solver.solve(&rhs).expect("solve should succeed");

        for i in 0..n {
            assert!(
                (x[i] - rhs[i]).abs() < 1e-14,
                "identity solve failed at index {i}: got {} expected {}",
                x[i],
                rhs[i]
            );
        }
    }

    #[test]
    fn test_nnz_factor_count() {
        let (row_offsets, col_indices, values, n) = spd_3x3_lower();

        let mut solver = SupernodalCholeskySolver::symbolic(&row_offsets, &col_indices, n)
            .expect("symbolic should succeed");
        solver
            .numeric(&row_offsets, &col_indices, &values)
            .expect("numeric should succeed");

        let nnz = solver.nnz_factor();
        // For 3x3 tridiagonal:
        // L has entries at (0,0), (1,0), (1,1), (2,1), (2,2) = 5 entries
        assert!(nnz >= 5, "nnz_factor = {nnz}, expected >= 5");
    }
}