scirs2-sparse 0.4.2

Sparse matrix module for SciRS2 (scirs2-sparse)
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
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
//! Matrix decomposition algorithms for sparse matrices
//!
//! This module provides various matrix decomposition algorithms optimized
//! for sparse matrices, including LU, QR, Cholesky, and incomplete variants.

use crate::csr_array::CsrArray;
use crate::error::{SparseError, SparseResult};
use crate::sparray::SparseArray;
use scirs2_core::ndarray::{Array1, Array2};
use scirs2_core::numeric::Float;
use scirs2_core::SparseElement;
use std::collections::HashMap;
use std::fmt::Debug;
use std::ops::{Add, Div, Mul, Sub};

/// LU decomposition result
#[derive(Debug, Clone)]
pub struct LUResult<T>
where
    T: Float + SparseElement + Debug + Copy + 'static,
{
    /// Lower triangular factor
    pub l: CsrArray<T>,
    /// Upper triangular factor
    pub u: CsrArray<T>,
    /// Permutation matrix (as permutation vector)
    pub p: Array1<usize>,
    /// Whether decomposition was successful
    pub success: bool,
}

/// QR decomposition result
#[derive(Debug, Clone)]
pub struct QRResult<T>
where
    T: Float + SparseElement + Debug + Copy + 'static,
{
    /// Orthogonal factor Q
    pub q: CsrArray<T>,
    /// Upper triangular factor R
    pub r: CsrArray<T>,
    /// Whether decomposition was successful
    pub success: bool,
}

/// Cholesky decomposition result
#[derive(Debug, Clone)]
pub struct CholeskyResult<T>
where
    T: Float + SparseElement + Debug + Copy + 'static,
{
    /// Lower triangular Cholesky factor
    pub l: CsrArray<T>,
    /// Whether decomposition was successful
    pub success: bool,
}

/// Pivoted Cholesky decomposition result
#[derive(Debug, Clone)]
pub struct PivotedCholeskyResult<T>
where
    T: Float + SparseElement + Debug + Copy + 'static,
{
    /// Lower triangular Cholesky factor
    pub l: CsrArray<T>,
    /// Permutation matrix (as permutation vector)
    pub p: Array1<usize>,
    /// Rank of the decomposition (number of positive eigenvalues)
    pub rank: usize,
    /// Whether decomposition was successful
    pub success: bool,
}

/// Pivoting strategy for LU decomposition
#[derive(Debug, Clone, Default)]
pub enum PivotingStrategy {
    /// No pivoting (fastest but potentially unstable)
    None,
    /// Partial pivoting - choose largest element in column (default)
    #[default]
    Partial,
    /// Threshold pivoting - partial pivoting with threshold
    Threshold(f64),
    /// Scaled partial pivoting - account for row scaling
    ScaledPartial,
    /// Complete pivoting - choose largest element in submatrix (most stable but expensive)
    Complete,
    /// Rook pivoting - hybrid approach balancing stability and cost
    Rook,
}

/// Options for LU decomposition
#[derive(Debug, Clone)]
pub struct LUOptions {
    /// Pivoting strategy to use
    pub pivoting: PivotingStrategy,
    /// Threshold for numerical zero (default: 1e-14)
    pub zero_threshold: f64,
    /// Whether to check for singularity (default: true)
    pub check_singular: bool,
}

impl Default for LUOptions {
    fn default() -> Self {
        Self {
            pivoting: PivotingStrategy::default(),
            zero_threshold: 1e-14,
            check_singular: true,
        }
    }
}

/// Options for incomplete LU decomposition
#[derive(Debug, Clone)]
pub struct ILUOptions {
    /// Drop tolerance for numerical stability
    pub drop_tol: f64,
    /// Fill factor (maximum fill-in ratio)
    pub fill_factor: f64,
    /// Maximum number of fill-in entries per row
    pub max_fill_per_row: usize,
    /// Pivoting strategy to use
    pub pivoting: PivotingStrategy,
}

impl Default for ILUOptions {
    fn default() -> Self {
        Self {
            drop_tol: 1e-4,
            fill_factor: 2.0,
            max_fill_per_row: 20,
            pivoting: PivotingStrategy::default(),
        }
    }
}

/// Options for incomplete Cholesky decomposition
#[derive(Debug, Clone)]
pub struct ICOptions {
    /// Drop tolerance for numerical stability
    pub drop_tol: f64,
    /// Fill factor (maximum fill-in ratio)
    pub fill_factor: f64,
    /// Maximum number of fill-in entries per row
    pub max_fill_per_row: usize,
}

impl Default for ICOptions {
    fn default() -> Self {
        Self {
            drop_tol: 1e-4,
            fill_factor: 2.0,
            max_fill_per_row: 20,
        }
    }
}

/// Compute sparse LU decomposition with partial pivoting (backward compatibility)
///
/// Computes the LU decomposition of a sparse matrix A such that P*A = L*U,
/// where P is a permutation matrix, L is lower triangular, and U is upper triangular.
///
/// # Arguments
///
/// * `matrix` - The sparse matrix to decompose
/// * `pivot_threshold` - Pivoting threshold for numerical stability (0.0 to 1.0)
///
/// # Returns
///
/// LU decomposition result
///
/// # Examples
///
/// ```
/// use scirs2_sparse::linalg::lu_decomposition;
/// use scirs2_sparse::csr_array::CsrArray;
///
/// // Create a sparse matrix
/// let rows = vec![0, 0, 1, 2];
/// let cols = vec![0, 1, 1, 2];
/// let data = vec![2.0, 1.0, 3.0, 4.0];
/// let matrix = CsrArray::from_triplets(&rows, &cols, &data, (3, 3), false).expect("Operation failed");
///
/// let lu_result = lu_decomposition(&matrix, 0.1).expect("Operation failed");
/// ```
#[allow(dead_code)]
pub fn lu_decomposition<T, S>(_matrix: &S, pivotthreshold: f64) -> SparseResult<LUResult<T>>
where
    T: Float
        + SparseElement
        + Debug
        + Copy
        + Add<Output = T>
        + Sub<Output = T>
        + Mul<Output = T>
        + Div<Output = T>,
    S: SparseArray<T>,
{
    // Use _threshold pivoting for backward compatibility
    let options = LUOptions {
        pivoting: PivotingStrategy::Threshold(pivotthreshold),
        zero_threshold: 1e-14,
        check_singular: true,
    };

    lu_decomposition_with_options(_matrix, Some(options))
}

/// Compute sparse LU decomposition with enhanced pivoting strategies
///
/// Computes the LU decomposition of a sparse matrix A such that P*A = L*U,
/// where P is a permutation matrix, L is lower triangular, and U is upper triangular.
/// This version supports multiple pivoting strategies for enhanced numerical stability.
///
/// # Arguments
///
/// * `matrix` - The sparse matrix to decompose
/// * `options` - LU decomposition options (pivoting strategy, thresholds, etc.)
///
/// # Returns
///
/// LU decomposition result
///
/// # Examples
///
/// ```
/// use scirs2_sparse::linalg::{lu_decomposition_with_options, LUOptions, PivotingStrategy};
/// use scirs2_sparse::csr_array::CsrArray;
///
/// // Create a sparse matrix
/// let rows = vec![0, 0, 1, 2];
/// let cols = vec![0, 1, 1, 2];
/// let data = vec![2.0, 1.0, 3.0, 4.0];
/// let matrix = CsrArray::from_triplets(&rows, &cols, &data, (3, 3), false).expect("Operation failed");
///
/// let options = LUOptions {
///     pivoting: PivotingStrategy::ScaledPartial,
///     zero_threshold: 1e-12,
///     check_singular: true,
/// };
/// let lu_result = lu_decomposition_with_options(&matrix, Some(options)).expect("Operation failed");
/// ```
#[allow(dead_code)]
pub fn lu_decomposition_with_options<T, S>(
    matrix: &S,
    options: Option<LUOptions>,
) -> SparseResult<LUResult<T>>
where
    T: Float
        + SparseElement
        + Debug
        + Copy
        + Add<Output = T>
        + Sub<Output = T>
        + Mul<Output = T>
        + Div<Output = T>,
    S: SparseArray<T>,
{
    let opts = options.unwrap_or_default();
    let (n, m) = matrix.shape();
    if n != m {
        return Err(SparseError::ValueError(
            "Matrix must be square for LU decomposition".to_string(),
        ));
    }

    // Convert to working format
    let (row_indices, col_indices, values) = matrix.find();
    let mut working_matrix = SparseWorkingMatrix::from_triplets(
        row_indices.as_slice().expect("Operation failed"),
        col_indices.as_slice().expect("Operation failed"),
        values.as_slice().expect("Operation failed"),
        n,
    );

    // Initialize permutations
    let mut row_perm: Vec<usize> = (0..n).collect();
    let mut col_perm: Vec<usize> = (0..n).collect();

    // Compute row scaling factors for scaled partial pivoting
    let mut row_scales = vec![T::sparse_one(); n];
    if matches!(opts.pivoting, PivotingStrategy::ScaledPartial) {
        for (i, scale) in row_scales.iter_mut().enumerate().take(n) {
            let row_data = working_matrix.get_row(i);
            let max_val = row_data
                .values()
                .map(|&v| v.abs())
                .fold(T::sparse_zero(), |a, b| if a > b { a } else { b });
            if max_val > T::sparse_zero() {
                *scale = max_val;
            }
        }
    }

    // Gaussian elimination with enhanced pivoting
    for k in 0..n - 1 {
        // Find pivot using selected strategy
        let (pivot_row, pivot_col) =
            find_enhanced_pivot(&working_matrix, k, &row_perm, &col_perm, &row_scales, &opts)?;

        // Apply row and column permutations
        if pivot_row != k {
            row_perm.swap(k, pivot_row);
        }
        if pivot_col != k
            && matches!(
                opts.pivoting,
                PivotingStrategy::Complete | PivotingStrategy::Rook
            )
        {
            col_perm.swap(k, pivot_col);
            // When columns are swapped, we need to update all matrix elements
            for &row_idx in row_perm.iter().take(n) {
                let temp = working_matrix.get(row_idx, k);
                working_matrix.set(row_idx, k, working_matrix.get(row_idx, pivot_col));
                working_matrix.set(row_idx, pivot_col, temp);
            }
        }

        let actual_pivot_row = row_perm[k];
        let actual_pivot_col = col_perm[k];
        let pivot_value = working_matrix.get(actual_pivot_row, actual_pivot_col);

        // Check for numerical singularity
        if opts.check_singular
            && pivot_value.abs() < T::from(opts.zero_threshold).expect("Operation failed")
        {
            return Ok(LUResult {
                l: CsrArray::from_triplets(&[], &[], &[], (n, n), false)?,
                u: CsrArray::from_triplets(&[], &[], &[], (n, n), false)?,
                p: Array1::from_vec(row_perm),
                success: false,
            });
        }

        // Eliminate below pivot
        for &actual_row_i in row_perm.iter().take(n).skip(k + 1) {
            let factor = working_matrix.get(actual_row_i, actual_pivot_col) / pivot_value;

            if !SparseElement::is_zero(&factor) {
                // Store multiplier in L
                working_matrix.set(actual_row_i, actual_pivot_col, factor);

                // Update row i
                let pivot_row_data = working_matrix.get_row(actual_pivot_row);
                for (col, &value) in &pivot_row_data {
                    if *col > k {
                        let old_val = working_matrix.get(actual_row_i, *col);
                        working_matrix.set(actual_row_i, *col, old_val - factor * value);
                    }
                }
            }
        }
    }

    // Extract L and U matrices with proper permutation
    let (l_rows, l_cols, l_vals, u_rows, u_cols, u_vals) =
        extract_lu_factors(&working_matrix, &row_perm, n);

    let l = CsrArray::from_triplets(&l_rows, &l_cols, &l_vals, (n, n), false)?;
    let u = CsrArray::from_triplets(&u_rows, &u_cols, &u_vals, (n, n), false)?;

    Ok(LUResult {
        l,
        u,
        p: Array1::from_vec(row_perm),
        success: true,
    })
}

/// Compute sparse QR decomposition using Givens rotations
///
/// Computes the QR decomposition of a sparse matrix A = Q*R,
/// where Q is orthogonal and R is upper triangular.
///
/// # Arguments
///
/// * `matrix` - The sparse matrix to decompose
///
/// # Returns
///
/// QR decomposition result
///
/// # Examples
///
/// ```
/// use scirs2_sparse::linalg::qr_decomposition;
/// use scirs2_sparse::csr_array::CsrArray;
///
/// let rows = vec![0, 1, 2];
/// let cols = vec![0, 0, 1];
/// let data = vec![1.0, 2.0, 3.0];
/// let matrix = CsrArray::from_triplets(&rows, &cols, &data, (3, 2), false).expect("Operation failed");
///
/// let qr_result = qr_decomposition(&matrix).expect("Operation failed");
/// ```
#[allow(dead_code)]
pub fn qr_decomposition<T, S>(matrix: &S) -> SparseResult<QRResult<T>>
where
    T: Float
        + SparseElement
        + Debug
        + Copy
        + Add<Output = T>
        + Sub<Output = T>
        + Mul<Output = T>
        + Div<Output = T>,
    S: SparseArray<T>,
{
    let (m, n) = matrix.shape();

    // Convert to dense for QR (sparse QR is complex)
    let dense_matrix = matrix.to_array();

    // Simple Gram-Schmidt QR decomposition
    let mut q = Array2::zeros((m, n));
    let mut r = Array2::zeros((n, n));

    for j in 0..n {
        // Copy column j
        for i in 0..m {
            q[[i, j]] = dense_matrix[[i, j]];
        }

        // Orthogonalize against previous columns
        for k in 0..j {
            let mut dot = T::sparse_zero();
            for i in 0..m {
                dot = dot + q[[i, k]] * dense_matrix[[i, j]];
            }
            r[[k, j]] = dot;

            for i in 0..m {
                q[[i, j]] = q[[i, j]] - dot * q[[i, k]];
            }
        }

        // Normalize
        let mut norm = T::sparse_zero();
        for i in 0..m {
            norm = norm + q[[i, j]] * q[[i, j]];
        }
        norm = norm.sqrt();
        r[[j, j]] = norm;

        if !SparseElement::is_zero(&norm) {
            for i in 0..m {
                q[[i, j]] = q[[i, j]] / norm;
            }
        }
    }

    // Convert back to sparse
    let q_sparse = dense_to_sparse(&q)?;
    let r_sparse = dense_to_sparse(&r)?;

    Ok(QRResult {
        q: q_sparse,
        r: r_sparse,
        success: true,
    })
}

/// Compute sparse Cholesky decomposition
///
/// Computes the Cholesky decomposition of a symmetric positive definite matrix A = L*L^T,
/// where L is lower triangular.
///
/// # Arguments
///
/// * `matrix` - The symmetric positive definite sparse matrix
///
/// # Returns
///
/// Cholesky decomposition result
///
/// # Examples
///
/// ```
/// use scirs2_sparse::linalg::cholesky_decomposition;
/// use scirs2_sparse::csr_array::CsrArray;
///
/// // Create a simple SPD matrix
/// let rows = vec![0, 1, 1, 2, 2, 2];
/// let cols = vec![0, 0, 1, 0, 1, 2];
/// let data = vec![4.0, 2.0, 5.0, 1.0, 3.0, 6.0];
/// let matrix = CsrArray::from_triplets(&rows, &cols, &data, (3, 3), false).expect("Operation failed");
///
/// let chol_result = cholesky_decomposition(&matrix).expect("Operation failed");
/// ```
#[allow(dead_code)]
pub fn cholesky_decomposition<T, S>(matrix: &S) -> SparseResult<CholeskyResult<T>>
where
    T: Float
        + SparseElement
        + Debug
        + Copy
        + Add<Output = T>
        + Sub<Output = T>
        + Mul<Output = T>
        + Div<Output = T>,
    S: SparseArray<T>,
{
    let (n, m) = matrix.shape();
    if n != m {
        return Err(SparseError::ValueError(
            "Matrix must be square for Cholesky decomposition".to_string(),
        ));
    }

    // Convert to working format
    let (row_indices, col_indices, values) = matrix.find();
    let mut working_matrix = SparseWorkingMatrix::from_triplets(
        row_indices.as_slice().expect("Operation failed"),
        col_indices.as_slice().expect("Operation failed"),
        values.as_slice().expect("Operation failed"),
        n,
    );

    // Cholesky decomposition algorithm
    for k in 0..n {
        // Compute diagonal element
        let mut sum = T::sparse_zero();
        for j in 0..k {
            let l_kj = working_matrix.get(k, j);
            sum = sum + l_kj * l_kj;
        }

        let a_kk = working_matrix.get(k, k);
        let diag_val = a_kk - sum;

        if diag_val <= T::sparse_zero() {
            return Ok(CholeskyResult {
                l: CsrArray::from_triplets(&[], &[], &[], (n, n), false)?,
                success: false,
            });
        }

        let l_kk = diag_val.sqrt();
        working_matrix.set(k, k, l_kk);

        // Compute below-diagonal elements
        for i in (k + 1)..n {
            let mut sum = T::sparse_zero();
            for j in 0..k {
                sum = sum + working_matrix.get(i, j) * working_matrix.get(k, j);
            }

            let a_ik = working_matrix.get(i, k);
            let l_ik = (a_ik - sum) / l_kk;
            working_matrix.set(i, k, l_ik);
        }
    }

    // Extract lower triangular _matrix
    let (l_rows, l_cols, l_vals) = extract_lower_triangular(&working_matrix, n);
    let l = CsrArray::from_triplets(&l_rows, &l_cols, &l_vals, (n, n), false)?;

    Ok(CholeskyResult { l, success: true })
}

/// Compute pivoted Cholesky decomposition
///
/// Computes the pivoted Cholesky decomposition of a symmetric matrix A = P^T * L * L^T * P,
/// where P is a permutation matrix and L is lower triangular. This version can handle
/// indefinite matrices by determining the rank and producing a partial decomposition.
///
/// # Arguments
///
/// * `matrix` - The symmetric sparse matrix
/// * `threshold` - Pivoting threshold for numerical stability (default: 1e-12)
///
/// # Returns
///
/// Pivoted Cholesky decomposition result with rank determination
///
/// # Examples
///
/// ```
/// use scirs2_sparse::linalg::pivoted_cholesky_decomposition;
/// use scirs2_sparse::csr_array::CsrArray;
///
/// // Create a symmetric indefinite matrix
/// let rows = vec![0, 1, 1, 2, 2, 2];
/// let cols = vec![0, 0, 1, 0, 1, 2];  
/// let data = vec![1.0, 2.0, -1.0, 3.0, 1.0, 2.0];
/// let matrix = CsrArray::from_triplets(&rows, &cols, &data, (3, 3), false).expect("Operation failed");
///
/// let chol_result = pivoted_cholesky_decomposition(&matrix, Some(1e-12)).expect("Operation failed");
/// ```
#[allow(dead_code)]
pub fn pivoted_cholesky_decomposition<T, S>(
    matrix: &S,
    threshold: Option<T>,
) -> SparseResult<PivotedCholeskyResult<T>>
where
    T: Float
        + SparseElement
        + Debug
        + Copy
        + Add<Output = T>
        + Sub<Output = T>
        + Mul<Output = T>
        + Div<Output = T>,
    S: SparseArray<T>,
{
    let (n, m) = matrix.shape();
    if n != m {
        return Err(SparseError::ValueError(
            "Matrix must be square for Cholesky decomposition".to_string(),
        ));
    }

    let threshold = threshold.unwrap_or_else(|| T::from(1e-12).expect("Operation failed"));

    // Convert to working format
    let (row_indices, col_indices, values) = matrix.find();
    let mut working_matrix = SparseWorkingMatrix::from_triplets(
        row_indices.as_slice().expect("Operation failed"),
        col_indices.as_slice().expect("Operation failed"),
        values.as_slice().expect("Operation failed"),
        n,
    );

    // Initialize permutation
    let mut perm: Vec<usize> = (0..n).collect();
    let mut rank = 0;

    // Pivoted Cholesky algorithm
    for k in 0..n {
        // Find the pivot: largest diagonal element among remaining
        let mut max_diag = T::sparse_zero();
        let mut pivot_idx = k;

        for i in k..n {
            let mut diag_val = working_matrix.get(perm[i], perm[i]);
            for j in 0..k {
                let l_ij = working_matrix.get(perm[i], perm[j]);
                diag_val = diag_val - l_ij * l_ij;
            }
            if diag_val > max_diag {
                max_diag = diag_val;
                pivot_idx = i;
            }
        }

        // Check if we should stop (matrix is not positive definite beyond this point)
        if max_diag <= threshold {
            break;
        }

        // Swap rows/columns in permutation
        if pivot_idx != k {
            perm.swap(k, pivot_idx);
        }

        // Compute L[k,k]
        let l_kk = max_diag.sqrt();
        working_matrix.set(perm[k], perm[k], l_kk);
        rank += 1;

        // Update column k below diagonal
        for i in (k + 1)..n {
            let mut sum = T::sparse_zero();
            for j in 0..k {
                sum = sum
                    + working_matrix.get(perm[i], perm[j]) * working_matrix.get(perm[k], perm[j]);
            }

            let a_ik = working_matrix.get(perm[i], perm[k]);
            let l_ik = (a_ik - sum) / l_kk;
            working_matrix.set(perm[i], perm[k], l_ik);
        }
    }

    // Extract lower triangular matrix with proper permutation
    let mut l_rows = Vec::new();
    let mut l_cols = Vec::new();
    let mut l_vals = Vec::new();

    for i in 0..rank {
        for j in 0..=i {
            let val = working_matrix.get(perm[i], perm[j]);
            if val != T::sparse_zero() {
                l_rows.push(i);
                l_cols.push(j);
                l_vals.push(val);
            }
        }
    }

    let l = CsrArray::from_triplets(&l_rows, &l_cols, &l_vals, (n, rank), false)?;
    let p = Array1::from_vec(perm);

    Ok(PivotedCholeskyResult {
        l,
        p,
        rank,
        success: true,
    })
}

/// LDLT decomposition result for symmetric indefinite matrices
#[derive(Debug, Clone)]
pub struct LDLTResult<T>
where
    T: Float + SparseElement + Debug + Copy + 'static,
{
    /// Lower triangular factor L (unit diagonal)
    pub l: CsrArray<T>,
    /// Diagonal factor D
    pub d: Array1<T>,
    /// Permutation matrix (as permutation vector)
    pub p: Array1<usize>,
    /// Whether decomposition was successful
    pub success: bool,
}

/// Compute LDLT decomposition for symmetric indefinite matrices
///
/// Computes the LDLT decomposition of a symmetric matrix A = P^T * L * D * L^T * P,
/// where P is a permutation matrix, L is unit lower triangular, and D is diagonal.
/// This method can handle indefinite matrices unlike Cholesky decomposition.
///
/// # Arguments
///
/// * `matrix` - The symmetric sparse matrix
/// * `pivoting` - Whether to use pivoting for numerical stability (default: true)
/// * `threshold` - Pivoting threshold for numerical stability (default: 1e-12)
///
/// # Returns
///
/// LDLT decomposition result
///
/// # Examples
///
/// ```
/// use scirs2_sparse::linalg::ldlt_decomposition;
/// use scirs2_sparse::csr_array::CsrArray;
///
/// // Create a symmetric indefinite matrix
/// let rows = vec![0, 1, 1, 2, 2, 2];
/// let cols = vec![0, 0, 1, 0, 1, 2];  
/// let data = vec![1.0, 2.0, -1.0, 3.0, 1.0, 2.0];
/// let matrix = CsrArray::from_triplets(&rows, &cols, &data, (3, 3), false).expect("Operation failed");
///
/// let ldlt_result = ldlt_decomposition(&matrix, Some(true), Some(1e-12)).expect("Operation failed");
/// ```
#[allow(dead_code)]
pub fn ldlt_decomposition<T, S>(
    matrix: &S,
    pivoting: Option<bool>,
    threshold: Option<T>,
) -> SparseResult<LDLTResult<T>>
where
    T: Float
        + SparseElement
        + Debug
        + Copy
        + Add<Output = T>
        + Sub<Output = T>
        + Mul<Output = T>
        + Div<Output = T>,
    S: SparseArray<T>,
{
    let (n, m) = matrix.shape();
    if n != m {
        return Err(SparseError::ValueError(
            "Matrix must be square for LDLT decomposition".to_string(),
        ));
    }

    let use_pivoting = pivoting.unwrap_or(true);
    let threshold = threshold.unwrap_or_else(|| T::from(1e-12).expect("Operation failed"));

    // Convert to working format
    let (row_indices, col_indices, values) = matrix.find();
    let mut working_matrix = SparseWorkingMatrix::from_triplets(
        row_indices.as_slice().expect("Operation failed"),
        col_indices.as_slice().expect("Operation failed"),
        values.as_slice().expect("Operation failed"),
        n,
    );

    // Initialize permutation
    let mut perm: Vec<usize> = (0..n).collect();
    let mut d_values = vec![T::sparse_zero(); n];

    // LDLT decomposition with optional pivoting
    for k in 0..n {
        // Find pivot if pivoting is enabled
        if use_pivoting {
            let pivot_idx = find_ldlt_pivot(&working_matrix, k, &perm, threshold);
            if pivot_idx != k {
                perm.swap(k, pivot_idx);
            }
        }

        let actual_k = perm[k];

        // Compute diagonal element D[k,k]
        let mut diag_val = working_matrix.get(actual_k, actual_k);
        for j in 0..k {
            let l_kj = working_matrix.get(actual_k, perm[j]);
            diag_val = diag_val - l_kj * l_kj * d_values[j];
        }

        d_values[k] = diag_val;

        // Check for numerical issues
        if diag_val.abs() < threshold {
            return Ok(LDLTResult {
                l: CsrArray::from_triplets(&[], &[], &[], (n, n), false)?,
                d: Array1::from_vec(d_values),
                p: Array1::from_vec(perm),
                success: false,
            });
        }

        // Compute column k of L below the diagonal
        for i in (k + 1)..n {
            let actual_i = perm[i];
            let mut l_ik = working_matrix.get(actual_i, actual_k);

            for j in 0..k {
                l_ik = l_ik
                    - working_matrix.get(actual_i, perm[j])
                        * working_matrix.get(actual_k, perm[j])
                        * d_values[j];
            }

            l_ik = l_ik / diag_val;
            working_matrix.set(actual_i, actual_k, l_ik);
        }

        // Set diagonal element of L to 1
        working_matrix.set(actual_k, actual_k, T::sparse_one());
    }

    // Extract L matrix (unit lower triangular)
    let (l_rows, l_cols, l_vals) = extract_unit_lower_triangular(&working_matrix, &perm, n);
    let l = CsrArray::from_triplets(&l_rows, &l_cols, &l_vals, (n, n), false)?;

    Ok(LDLTResult {
        l,
        d: Array1::from_vec(d_values),
        p: Array1::from_vec(perm),
        success: true,
    })
}

/// Find pivot for LDLT decomposition using Bunch-Kaufman strategy
#[allow(dead_code)]
fn find_ldlt_pivot<T>(
    matrix: &SparseWorkingMatrix<T>,
    k: usize,
    perm: &[usize],
    threshold: T,
) -> usize
where
    T: Float + SparseElement + Debug + Copy,
{
    let n = matrix.n;
    let mut max_val = T::sparse_zero();
    let mut pivot_idx = k;

    // Look for largest diagonal element among remaining rows
    for (i, &actual_i) in perm.iter().enumerate().take(n).skip(k) {
        let diag_val = matrix.get(actual_i, actual_i).abs();

        if diag_val > max_val {
            max_val = diag_val;
            pivot_idx = i;
        }
    }

    // Check if pivot is acceptable
    if max_val >= threshold {
        pivot_idx
    } else {
        k // Use current position if no good pivot found
    }
}

/// Extract unit lower triangular matrix from working matrix
#[allow(dead_code)]
fn extract_unit_lower_triangular<T>(
    matrix: &SparseWorkingMatrix<T>,
    perm: &[usize],
    n: usize,
) -> (Vec<usize>, Vec<usize>, Vec<T>)
where
    T: Float + SparseElement + Debug + Copy,
{
    let mut rows = Vec::new();
    let mut cols = Vec::new();
    let mut vals = Vec::new();

    for i in 0..n {
        let actual_i = perm[i];

        // Add diagonal element (always 1 for unit triangular)
        rows.push(i);
        cols.push(i);
        vals.push(T::sparse_one());

        // Add below-diagonal elements
        for (j, &perm_j) in perm.iter().enumerate().take(i) {
            let val = matrix.get(actual_i, perm_j);
            if val != T::sparse_zero() {
                rows.push(i);
                cols.push(j);
                vals.push(val);
            }
        }
    }

    (rows, cols, vals)
}

/// Compute incomplete LU decomposition (ILU)
///
/// Computes an approximate LU decomposition with controlled fill-in
/// for use as a preconditioner in iterative methods.
///
/// # Arguments
///
/// * `matrix` - The sparse matrix to decompose
/// * `options` - ILU options controlling fill-in and dropping
///
/// # Returns
///
/// Incomplete LU decomposition result
#[allow(dead_code)]
pub fn incomplete_lu<T, S>(matrix: &S, options: Option<ILUOptions>) -> SparseResult<LUResult<T>>
where
    T: Float
        + SparseElement
        + Debug
        + Copy
        + Add<Output = T>
        + Sub<Output = T>
        + Mul<Output = T>
        + Div<Output = T>,
    S: SparseArray<T>,
{
    let opts = options.unwrap_or_default();
    let (n, m) = matrix.shape();

    if n != m {
        return Err(SparseError::ValueError(
            "Matrix must be square for ILU decomposition".to_string(),
        ));
    }

    // Convert to working format
    let (row_indices, col_indices, values) = matrix.find();
    let mut working_matrix = SparseWorkingMatrix::from_triplets(
        row_indices.as_slice().expect("Operation failed"),
        col_indices.as_slice().expect("Operation failed"),
        values.as_slice().expect("Operation failed"),
        n,
    );

    // ILU(0) algorithm - no fill-in beyond original sparsity pattern
    for k in 0..n - 1 {
        let pivot_val = working_matrix.get(k, k);

        if pivot_val.abs() < T::from(1e-14).expect("Operation failed") {
            continue; // Skip singular pivot
        }

        // Get all non-zero entries in column k below diagonal
        let col_k_entries = working_matrix.get_column_below_diagonal(k);

        for &row_i in &col_k_entries {
            let factor = working_matrix.get(row_i, k) / pivot_val;

            // Drop small factors
            if factor.abs() < T::from(opts.drop_tol).expect("Operation failed") {
                working_matrix.set(row_i, k, T::sparse_zero());
                continue;
            }

            working_matrix.set(row_i, k, factor);

            // Update row i (only existing non-zeros)
            let row_k_entries = working_matrix.get_row_after_column(k, k);
            for (col_j, &val_kj) in &row_k_entries {
                if working_matrix.has_entry(row_i, *col_j) {
                    let old_val = working_matrix.get(row_i, *col_j);
                    let new_val = old_val - factor * val_kj;

                    // Drop small values
                    if new_val.abs() < T::from(opts.drop_tol).expect("Operation failed") {
                        working_matrix.set(row_i, *col_j, T::sparse_zero());
                    } else {
                        working_matrix.set(row_i, *col_j, new_val);
                    }
                }
            }
        }
    }

    // Extract L and U factors
    let identity_p: Vec<usize> = (0..n).collect();
    let (l_rows, l_cols, l_vals, u_rows, u_cols, u_vals) =
        extract_lu_factors(&working_matrix, &identity_p, n);

    let l = CsrArray::from_triplets(&l_rows, &l_cols, &l_vals, (n, n), false)?;
    let u = CsrArray::from_triplets(&u_rows, &u_cols, &u_vals, (n, n), false)?;

    Ok(LUResult {
        l,
        u,
        p: Array1::from_vec(identity_p),
        success: true,
    })
}

/// Compute incomplete Cholesky decomposition (IC)
///
/// Computes an approximate Cholesky decomposition with controlled fill-in
/// for use as a preconditioner in iterative methods.
///
/// # Arguments
///
/// * `matrix` - The symmetric positive definite sparse matrix
/// * `options` - IC options controlling fill-in and dropping
///
/// # Returns
///
/// Incomplete Cholesky decomposition result
#[allow(dead_code)]
pub fn incomplete_cholesky<T, S>(
    matrix: &S,
    options: Option<ICOptions>,
) -> SparseResult<CholeskyResult<T>>
where
    T: Float
        + SparseElement
        + Debug
        + Copy
        + Add<Output = T>
        + Sub<Output = T>
        + Mul<Output = T>
        + Div<Output = T>,
    S: SparseArray<T>,
{
    let opts = options.unwrap_or_default();
    let (n, m) = matrix.shape();

    if n != m {
        return Err(SparseError::ValueError(
            "Matrix must be square for IC decomposition".to_string(),
        ));
    }

    // Convert to working format
    let (row_indices, col_indices, values) = matrix.find();
    let mut working_matrix = SparseWorkingMatrix::from_triplets(
        row_indices.as_slice().expect("Operation failed"),
        col_indices.as_slice().expect("Operation failed"),
        values.as_slice().expect("Operation failed"),
        n,
    );

    // IC(0) algorithm - no fill-in beyond original sparsity pattern
    for k in 0..n {
        // Compute diagonal element
        let mut sum = T::sparse_zero();
        let row_k_before_k = working_matrix.get_row_before_column(k, k);
        for &val_kj in row_k_before_k.values() {
            sum = sum + val_kj * val_kj;
        }

        let a_kk = working_matrix.get(k, k);
        let diag_val = a_kk - sum;

        if diag_val <= T::sparse_zero() {
            return Ok(CholeskyResult {
                l: CsrArray::from_triplets(&[], &[], &[], (n, n), false)?,
                success: false,
            });
        }

        let l_kk = diag_val.sqrt();
        working_matrix.set(k, k, l_kk);

        // Compute below-diagonal elements (only existing entries)
        let col_k_below = working_matrix.get_column_below_diagonal(k);
        for &row_i in &col_k_below {
            let mut sum = T::sparse_zero();
            let row_i_before_k = working_matrix.get_row_before_column(row_i, k);
            let row_k_before_k = working_matrix.get_row_before_column(k, k);

            // Compute dot product of L[i, :k] and L[k, :k]
            for (col_j, &val_ij) in &row_i_before_k {
                if let Some(&val_kj) = row_k_before_k.get(col_j) {
                    sum = sum + val_ij * val_kj;
                }
            }

            let a_ik = working_matrix.get(row_i, k);
            let l_ik = (a_ik - sum) / l_kk;

            // Drop small values
            if l_ik.abs() < T::from(opts.drop_tol).expect("Operation failed") {
                working_matrix.set(row_i, k, T::sparse_zero());
            } else {
                working_matrix.set(row_i, k, l_ik);
            }
        }
    }

    // Extract lower triangular matrix
    let (l_rows, l_cols, l_vals) = extract_lower_triangular(&working_matrix, n);
    let l = CsrArray::from_triplets(&l_rows, &l_cols, &l_vals, (n, n), false)?;

    Ok(CholeskyResult { l, success: true })
}

/// Simple sparse working matrix for decomposition algorithms
struct SparseWorkingMatrix<T>
where
    T: Float + SparseElement + Debug + Copy,
{
    data: HashMap<(usize, usize), T>,
    n: usize,
}

impl<T> SparseWorkingMatrix<T>
where
    T: Float
        + SparseElement
        + Debug
        + Copy
        + Add<Output = T>
        + Sub<Output = T>
        + Mul<Output = T>
        + Div<Output = T>,
{
    fn from_triplets(rows: &[usize], cols: &[usize], values: &[T], n: usize) -> Self {
        let mut data = HashMap::new();

        for (i, (&row, &col)) in rows.iter().zip(cols.iter()).enumerate() {
            data.insert((row, col), values[i]);
        }

        Self { data, n }
    }

    fn get(&self, row: usize, col: usize) -> T {
        self.data
            .get(&(row, col))
            .copied()
            .unwrap_or(T::sparse_zero())
    }

    fn set(&mut self, row: usize, col: usize, value: T) {
        if SparseElement::is_zero(&value) {
            self.data.remove(&(row, col));
        } else {
            self.data.insert((row, col), value);
        }
    }

    fn has_entry(&self, row: usize, col: usize) -> bool {
        self.data.contains_key(&(row, col))
    }

    fn get_row(&self, row: usize) -> HashMap<usize, T> {
        let mut result = HashMap::new();
        for (&(r, c), &value) in &self.data {
            if r == row {
                result.insert(c, value);
            }
        }
        result
    }

    fn get_row_after_column(&self, row: usize, col: usize) -> HashMap<usize, T> {
        let mut result = HashMap::new();
        for (&(r, c), &value) in &self.data {
            if r == row && c > col {
                result.insert(c, value);
            }
        }
        result
    }

    fn get_row_before_column(&self, row: usize, col: usize) -> HashMap<usize, T> {
        let mut result = HashMap::new();
        for (&(r, c), &value) in &self.data {
            if r == row && c < col {
                result.insert(c, value);
            }
        }
        result
    }

    fn get_column_below_diagonal(&self, col: usize) -> Vec<usize> {
        let mut result = Vec::new();
        for &(r, c) in self.data.keys() {
            if c == col && r > col {
                result.push(r);
            }
        }
        result.sort();
        result
    }
}

/// Find pivot for LU decomposition (backward compatibility)
#[allow(dead_code)]
fn find_pivot<T>(
    matrix: &SparseWorkingMatrix<T>,
    k: usize,
    p: &[usize],
    threshold: f64,
) -> SparseResult<usize>
where
    T: Float + SparseElement + Debug + Copy,
{
    // Use threshold pivoting for backward compatibility
    let opts = LUOptions {
        pivoting: PivotingStrategy::Threshold(threshold),
        zero_threshold: 1e-14,
        check_singular: true,
    };

    let row_scales = vec![T::sparse_one(); matrix.n];
    let col_perm: Vec<usize> = (0..matrix.n).collect();

    let (pivot_row, pivot_col) = find_enhanced_pivot(matrix, k, p, &col_perm, &row_scales, &opts)?;
    Ok(pivot_row)
}

/// Enhanced pivoting function supporting multiple strategies
#[allow(dead_code)]
fn find_enhanced_pivot<T>(
    matrix: &SparseWorkingMatrix<T>,
    k: usize,
    row_perm: &[usize],
    col_perm: &[usize],
    row_scales: &[T],
    opts: &LUOptions,
) -> SparseResult<(usize, usize)>
where
    T: Float + SparseElement + Debug + Copy,
{
    let n = matrix.n;

    match &opts.pivoting {
        PivotingStrategy::None => {
            // No pivoting - use diagonal element
            Ok((k, k))
        }

        PivotingStrategy::Partial => {
            // Standard partial pivoting - find largest element in column k
            let mut max_val = T::sparse_zero();
            let mut pivot_row = k;

            for (idx, &actual_row) in row_perm.iter().enumerate().skip(k).take(n - k) {
                let i = k + idx;
                let val = matrix.get(actual_row, col_perm[k]).abs();
                if val > max_val {
                    max_val = val;
                    pivot_row = i;
                }
            }

            Ok((pivot_row, k))
        }

        PivotingStrategy::Threshold(threshold) => {
            // Threshold pivoting - use first element above threshold
            let threshold_val = T::from(*threshold).expect("Operation failed");
            let mut max_val = T::sparse_zero();
            let mut pivot_row = k;

            for (idx, &actual_row) in row_perm.iter().enumerate().skip(k).take(n - k) {
                let i = k + idx;
                let val = matrix.get(actual_row, col_perm[k]).abs();
                if val > max_val {
                    max_val = val;
                    pivot_row = i;
                }
                // Use first element above threshold for efficiency
                if val >= threshold_val {
                    pivot_row = i;
                    break;
                }
            }

            Ok((pivot_row, k))
        }

        PivotingStrategy::ScaledPartial => {
            // Scaled partial pivoting - account for row scaling
            let mut max_ratio = T::sparse_zero();
            let mut pivot_row = k;

            for (idx, &actual_row) in row_perm.iter().enumerate().skip(k).take(n - k) {
                let i = k + idx;
                let val = matrix.get(actual_row, col_perm[k]).abs();
                let scale = row_scales[actual_row];

                let ratio = if scale > T::sparse_zero() {
                    val / scale
                } else {
                    val
                };

                if ratio > max_ratio {
                    max_ratio = ratio;
                    pivot_row = i;
                }
            }

            Ok((pivot_row, k))
        }

        PivotingStrategy::Complete => {
            // Complete pivoting - find largest element in remaining submatrix
            let mut max_val = T::sparse_zero();
            let mut pivot_row = k;
            let mut pivot_col = k;

            for (i_idx, &actual_row) in row_perm.iter().enumerate().skip(k).take(n - k) {
                let i = k + i_idx;
                for (j_idx, &actual_col) in col_perm.iter().enumerate().skip(k).take(n - k) {
                    let j = k + j_idx;
                    let val = matrix.get(actual_row, actual_col).abs();
                    if val > max_val {
                        max_val = val;
                        pivot_row = i;
                        pivot_col = j;
                    }
                }
            }

            Ok((pivot_row, pivot_col))
        }

        PivotingStrategy::Rook => {
            // Rook pivoting - alternating row and column searches
            let mut best_row = k;
            let mut best_col = k;
            let mut max_val = T::sparse_zero();

            // Start with partial pivoting in column k
            for (idx, &actual_row) in row_perm.iter().enumerate().skip(k).take(n - k) {
                let i = k + idx;
                let val = matrix.get(actual_row, col_perm[k]).abs();
                if val > max_val {
                    max_val = val;
                    best_row = i;
                }
            }

            // If we found a good pivot, check if we can improve by column pivoting
            if max_val > T::from(opts.zero_threshold).expect("Operation failed") {
                let actual_best_row = row_perm[best_row];
                let mut col_max = T::sparse_zero();

                for (idx, &actual_col) in col_perm.iter().enumerate().skip(k).take(n - k) {
                    let j = k + idx;
                    let val = matrix.get(actual_best_row, actual_col).abs();
                    if val > col_max {
                        col_max = val;
                        best_col = j;
                    }
                }

                // Use column pivot if it's significantly better
                let improvement_threshold = T::from(1.5).expect("Operation failed");
                if col_max > max_val * improvement_threshold {
                    // Recompute row pivot for the new column
                    max_val = T::sparse_zero();
                    for (idx, &actual_row) in row_perm.iter().enumerate().skip(k).take(n - k) {
                        let i = k + idx;
                        let val = matrix.get(actual_row, col_perm[best_col]).abs();
                        if val > max_val {
                            max_val = val;
                            best_row = i;
                        }
                    }
                }
            }

            Ok((best_row, best_col))
        }
    }
}

/// Extract L and U factors from working matrix
type LuFactors<T> = (
    Vec<usize>, // L row pointers
    Vec<usize>, // L column indices
    Vec<T>,     // L values
    Vec<usize>, // U row pointers
    Vec<usize>, // U column indices
    Vec<T>,     // U values
);

#[allow(dead_code)]
fn extract_lu_factors<T>(matrix: &SparseWorkingMatrix<T>, p: &[usize], n: usize) -> LuFactors<T>
where
    T: Float + SparseElement + Debug + Copy,
{
    let mut l_rows = Vec::new();
    let mut l_cols = Vec::new();
    let mut l_vals = Vec::new();
    let mut u_rows = Vec::new();
    let mut u_cols = Vec::new();
    let mut u_vals = Vec::new();

    #[allow(clippy::needless_range_loop)]
    for i in 0..n {
        let actual_row = p[i];

        // Add diagonal 1 to L
        l_rows.push(i);
        l_cols.push(i);
        l_vals.push(T::sparse_one());

        for j in 0..n {
            let val = matrix.get(actual_row, j);
            if !SparseElement::is_zero(&val) {
                if j < i {
                    // Below diagonal - goes to L
                    l_rows.push(i);
                    l_cols.push(j);
                    l_vals.push(val);
                } else {
                    // On or above diagonal - goes to U
                    u_rows.push(i);
                    u_cols.push(j);
                    u_vals.push(val);
                }
            }
        }
    }

    (l_rows, l_cols, l_vals, u_rows, u_cols, u_vals)
}

/// Extract lower triangular matrix
#[allow(dead_code)]
fn extract_lower_triangular<T>(
    matrix: &SparseWorkingMatrix<T>,
    n: usize,
) -> (Vec<usize>, Vec<usize>, Vec<T>)
where
    T: Float + SparseElement + Debug + Copy,
{
    let mut rows = Vec::new();
    let mut cols = Vec::new();
    let mut vals = Vec::new();

    for i in 0..n {
        for j in 0..=i {
            let val = matrix.get(i, j);
            if !SparseElement::is_zero(&val) {
                rows.push(i);
                cols.push(j);
                vals.push(val);
            }
        }
    }

    (rows, cols, vals)
}

/// Convert dense matrix to sparse
#[allow(dead_code)]
fn dense_to_sparse<T>(matrix: &Array2<T>) -> SparseResult<CsrArray<T>>
where
    T: Float + SparseElement + Debug + Copy,
{
    let (m, n) = matrix.dim();
    let mut rows = Vec::new();
    let mut cols = Vec::new();
    let mut vals = Vec::new();

    for i in 0..m {
        for j in 0..n {
            let val = matrix[[i, j]];
            if !SparseElement::is_zero(&val) {
                rows.push(i);
                cols.push(j);
                vals.push(val);
            }
        }
    }

    CsrArray::from_triplets(&rows, &cols, &vals, (m, n), false)
}

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

    fn create_test_matrix() -> CsrArray<f64> {
        // Create a simple test matrix
        let rows = vec![0, 0, 1, 1, 2, 2];
        let cols = vec![0, 1, 0, 1, 1, 2];
        let data = vec![2.0, 1.0, 1.0, 3.0, 2.0, 4.0];

        CsrArray::from_triplets(&rows, &cols, &data, (3, 3), false).expect("Operation failed")
    }

    fn create_spd_matrix() -> CsrArray<f64> {
        // Create a symmetric positive definite matrix
        let rows = vec![0, 1, 1, 2, 2, 2];
        let cols = vec![0, 0, 1, 0, 1, 2];
        let data = vec![4.0, 2.0, 5.0, 1.0, 3.0, 6.0];

        CsrArray::from_triplets(&rows, &cols, &data, (3, 3), false).expect("Operation failed")
    }

    #[test]
    fn test_lu_decomposition() {
        let matrix = create_test_matrix();
        let lu_result = lu_decomposition(&matrix, 0.1).expect("Operation failed");

        assert!(lu_result.success);
        assert_eq!(lu_result.l.shape(), (3, 3));
        assert_eq!(lu_result.u.shape(), (3, 3));
        assert_eq!(lu_result.p.len(), 3);
    }

    #[test]
    fn test_qr_decomposition() {
        let rows = vec![0, 1, 2];
        let cols = vec![0, 0, 1];
        let data = vec![1.0, 2.0, 3.0];
        let matrix =
            CsrArray::from_triplets(&rows, &cols, &data, (3, 2), false).expect("Operation failed");

        let qr_result = qr_decomposition(&matrix).expect("Operation failed");

        assert!(qr_result.success);
        assert_eq!(qr_result.q.shape(), (3, 2));
        assert_eq!(qr_result.r.shape(), (2, 2));
    }

    #[test]
    fn test_cholesky_decomposition() {
        let matrix = create_spd_matrix();
        let chol_result = cholesky_decomposition(&matrix).expect("Operation failed");

        assert!(chol_result.success);
        assert_eq!(chol_result.l.shape(), (3, 3));
    }

    #[test]
    fn test_incomplete_lu() {
        let matrix = create_test_matrix();
        let options = ILUOptions {
            drop_tol: 1e-6,
            ..Default::default()
        };

        let ilu_result = incomplete_lu(&matrix, Some(options)).expect("Operation failed");

        assert!(ilu_result.success);
        assert_eq!(ilu_result.l.shape(), (3, 3));
        assert_eq!(ilu_result.u.shape(), (3, 3));
    }

    #[test]
    fn test_incomplete_cholesky() {
        let matrix = create_spd_matrix();
        let options = ICOptions {
            drop_tol: 1e-6,
            ..Default::default()
        };

        let ic_result = incomplete_cholesky(&matrix, Some(options)).expect("Operation failed");

        assert!(ic_result.success);
        assert_eq!(ic_result.l.shape(), (3, 3));
    }

    #[test]
    fn test_sparse_working_matrix() {
        let rows = vec![0, 1, 2];
        let cols = vec![0, 1, 2];
        let vals = vec![1.0, 2.0, 3.0];

        let mut matrix = SparseWorkingMatrix::from_triplets(&rows, &cols, &vals, 3);

        assert_eq!(matrix.get(0, 0), 1.0);
        assert_eq!(matrix.get(1, 1), 2.0);
        assert_eq!(matrix.get(2, 2), 3.0);
        assert_eq!(matrix.get(0, 1), 0.0);

        matrix.set(0, 1, 5.0);
        assert_eq!(matrix.get(0, 1), 5.0);

        matrix.set(0, 1, 0.0);
        assert_eq!(matrix.get(0, 1), 0.0);
        assert!(!matrix.has_entry(0, 1));
    }

    #[test]
    fn test_dense_to_sparse_conversion() {
        let dense =
            Array2::from_shape_vec((2, 2), vec![1.0, 0.0, 2.0, 3.0]).expect("Operation failed");
        let sparse = dense_to_sparse(&dense).expect("Operation failed");

        assert_eq!(sparse.nnz(), 3);
        assert_eq!(sparse.get(0, 0), 1.0);
        assert_eq!(sparse.get(0, 1), 0.0);
        assert_eq!(sparse.get(1, 0), 2.0);
        assert_eq!(sparse.get(1, 1), 3.0);
    }
}