Skip to main content

scirs2_sparse/
lib.rs

1#![allow(clippy::manual_div_ceil)]
2#![allow(clippy::needless_return)]
3#![allow(clippy::manual_ok_err)]
4#![allow(clippy::needless_range_loop)]
5#![allow(clippy::while_let_loop)]
6#![allow(clippy::vec_init_then_push)]
7#![allow(clippy::should_implement_trait)]
8#![allow(clippy::only_used_in_recursion)]
9#![allow(clippy::manual_slice_fill)]
10#![allow(dead_code)]
11//! # SciRS2 Sparse - Sparse Matrix Operations
12//!
13//! **scirs2-sparse** provides comprehensive sparse matrix formats and operations modeled after SciPy's
14//! `sparse` module, offering CSR, CSC, COO, DOK, LIL, DIA, BSR formats with efficient algorithms
15//! for large-scale sparse linear algebra, eigenvalue problems, and graph operations.
16//!
17//! ## 🎯 Key Features
18//!
19//! - **SciPy Compatibility**: Drop-in replacement for `scipy.sparse` classes
20//! - **Multiple Formats**: CSR, CSC, COO, DOK, LIL, DIA, BSR with easy conversion
21//! - **Efficient Operations**: Sparse matrix-vector/matrix multiplication
22//! - **Linear Solvers**: Direct (LU, Cholesky) and iterative (CG, GMRES) solvers
23//! - **Eigenvalue Solvers**: ARPACK-based sparse eigenvalue computation
24//! - **Array API**: Modern NumPy-compatible array interface (recommended)
25//!
26//! ## 📦 Module Overview
27//!
28//! | SciRS2 Format | SciPy Equivalent | Description |
29//! |---------------|------------------|-------------|
30//! | `CsrArray` | `scipy.sparse.csr_array` | Compressed Sparse Row (efficient row slicing) |
31//! | `CscArray` | `scipy.sparse.csc_array` | Compressed Sparse Column (efficient column slicing) |
32//! | `CooArray` | `scipy.sparse.coo_array` | Coordinate format (efficient construction) |
33//! | `DokArray` | `scipy.sparse.dok_array` | Dictionary of Keys (efficient element access) |
34//! | `LilArray` | `scipy.sparse.lil_array` | List of Lists (efficient incremental construction) |
35//! | `DiaArray` | `scipy.sparse.dia_array` | Diagonal format (efficient banded matrices) |
36//! | `BsrArray` | `scipy.sparse.bsr_array` | Block Sparse Row (efficient block operations) |
37//!
38//! ## 🚀 Quick Start
39//!
40//! ```toml
41//! [dependencies]
42//! scirs2-sparse = "0.2.0"
43//! ```
44//!
45//! ```rust
46//! use scirs2_sparse::csr_array::CsrArray;
47//!
48//! // Create sparse matrix from triplets (row, col, value)
49//! let rows = vec![0, 0, 1, 2, 2];
50//! let cols = vec![0, 2, 2, 0, 1];
51//! let data = vec![1.0, 2.0, 3.0, 4.0, 5.0];
52//! let sparse = CsrArray::from_triplets(&rows, &cols, &data, (3, 3), false).expect("Operation failed");
53//! ```
54//!
55//! ## 🔒 Version: 0.2.0 (February 8, 2026)
56//!
57//! ## Matrix vs. Array API
58//!
59//! This module provides both a matrix-based API and an array-based API,
60//! following SciPy's transition to a more NumPy-compatible array interface.
61//!
62//! When using the array interface (e.g., `CsrArray`), please note that:
63//!
64//! - `*` performs element-wise multiplication, not matrix multiplication
65//! - Use `dot()` method for matrix multiplication
66//! - Operations like `sum` produce arrays, not matrices
67//! - Array-style slicing operations return scalars, 1D, or 2D arrays
68//!
69//! For new code, we recommend using the array interface, which is more consistent
70//! with the rest of the numerical ecosystem.
71//!
72//! ## Examples
73//!
74//! ### Matrix API (Legacy)
75//!
76//! ```
77//! use scirs2_sparse::csr::CsrMatrix;
78//!
79//! // Create a sparse matrix in CSR format
80//! let rows = vec![0, 0, 1, 2, 2];
81//! let cols = vec![0, 2, 2, 0, 1];
82//! let data = vec![1.0, 2.0, 3.0, 4.0, 5.0];
83//! let shape = (3, 3);
84//!
85//! let matrix = CsrMatrix::new(data, rows, cols, shape).expect("Operation failed");
86//! ```
87//!
88//! ### Array API (Recommended)
89//!
90//! ```
91//! use scirs2_sparse::csr_array::CsrArray;
92//!
93//! // Create a sparse array in CSR format
94//! let rows = vec![0, 0, 1, 2, 2];
95//! let cols = vec![0, 2, 2, 0, 1];
96//! let data = vec![1.0, 2.0, 3.0, 4.0, 5.0];
97//! let shape = (3, 3);
98//!
99//! // From triplets (COO-like construction)
100//! let array = CsrArray::from_triplets(&rows, &cols, &data, shape, false).expect("Operation failed");
101//!
102//! // Or directly from CSR components
103//! // let array = CsrArray::new(...);
104//! ```
105
106// Export error types
107pub mod error;
108pub use error::{SparseError, SparseResult};
109
110// Base trait for sparse arrays
111pub mod sparray;
112pub use sparray::{is_sparse, SparseArray, SparseSum};
113
114// Trait for symmetric sparse arrays
115pub mod sym_sparray;
116pub use sym_sparray::SymSparseArray;
117
118// No spatial module in sparse
119
120// Array API (recommended)
121pub mod csr_array;
122pub use csr_array::CsrArray;
123
124pub mod csc_array;
125pub use csc_array::CscArray;
126
127pub mod coo_array;
128pub use coo_array::CooArray;
129
130pub mod dok_array;
131pub use dok_array::DokArray;
132
133pub mod lil_array;
134pub use lil_array::LilArray;
135
136pub mod dia_array;
137pub use dia_array::DiaArray;
138
139pub mod bsr_array;
140pub use bsr_array::BsrArray;
141
142pub mod banded_array;
143pub use banded_array::BandedArray;
144
145// Symmetric array formats
146pub mod sym_csr;
147pub use sym_csr::{SymCsrArray, SymCsrMatrix};
148
149pub mod sym_coo;
150pub use sym_coo::{SymCooArray, SymCooMatrix};
151
152// Legacy matrix formats
153pub mod csr;
154pub use csr::CsrMatrix;
155
156pub mod csc;
157pub use csc::CscMatrix;
158
159pub mod coo;
160pub use coo::CooMatrix;
161
162pub mod dok;
163pub use dok::DokMatrix;
164
165pub mod lil;
166pub use lil::LilMatrix;
167
168pub mod dia;
169pub use dia::DiaMatrix;
170
171pub mod bsr;
172pub use bsr::BsrMatrix;
173
174pub mod banded;
175pub use banded::BandedMatrix;
176
177// Utility functions
178pub mod utils;
179
180// Linear algebra with sparse matrices
181pub mod linalg;
182// Re-export the main functions from the reorganized linalg module
183pub use linalg::{
184    // Functions from solvers
185    add,
186    // Functions from iterative
187    bicg,
188    bicgstab,
189    cg,
190    cholesky_decomposition,
191    // Enhanced operators
192    convolution_operator,
193    diag_matrix,
194    eigs,
195    eigsh,
196    enhanced_add,
197    enhanced_diagonal,
198    enhanced_scale,
199    enhanced_subtract,
200    expm,
201    // Functions from matfuncs
202    expm_multiply,
203    eye,
204    finite_difference_operator,
205    // GCROT solver
206    gcrot,
207    gmres,
208    incomplete_cholesky,
209    incomplete_lu,
210    inv,
211    // IRAM eigenvalue solver (v0.3.0)
212    iram,
213    iram_shift_invert,
214    lanczos,
215    // Decomposition functions
216    lu_decomposition,
217    matmul,
218    matrix_power,
219    multiply,
220    norm,
221    onenormest,
222    // Eigenvalue functions
223    power_iteration,
224    qr_decomposition,
225    // Specialized solvers (v0.2.0)
226    solve_arrow_matrix,
227    solve_banded_system,
228    solve_block_2x2,
229    solve_kronecker_system,
230    solve_saddle_point,
231    sparse_direct_solve,
232    sparse_lstsq,
233    spsolve,
234    svd_truncated,
235    // SVD functions
236    svds,
237    // TFQMR solver
238    tfqmr,
239    // IRAM (Arnoldi) eigenvalue solver (v0.3.0)
240    ArnoldiConfig,
241    ArpackOptions,
242    // Interfaces
243    AsLinearOperator,
244    // Types from iterative
245    BiCGOptions,
246    BiCGSTABOptions,
247    BiCGSTABResult,
248    // Enhanced operator types
249    BoundaryCondition,
250    CGOptions,
251    CGSOptions,
252    CGSResult,
253    CholeskyResult,
254    ConvolutionMode,
255    ConvolutionOperator,
256    // Operator types
257    DiagonalOperator,
258    EigenResult,
259    EigenvalueMethod,
260    EnhancedDiagonalOperator,
261    EnhancedDifferenceOperator,
262    EnhancedOperatorOptions,
263    EnhancedScaledOperator,
264    EnhancedSumOperator,
265    FiniteDifferenceOperator,
266    GCROTOptions,
267    GCROTResult,
268    GMRESOptions,
269    ICOptions,
270    // Preconditioners
271    ILU0Preconditioner,
272    ILUOptions,
273    IdentityOperator,
274    IterationResult,
275    JacobiPreconditioner,
276    // Decomposition types
277    LUResult,
278    LanczosOptions,
279    LinearOperator,
280    // Eigenvalue types
281    PowerIterationOptions,
282    QRResult,
283    SSORPreconditioner,
284    // SVD types
285    SVDOptions,
286    SVDResult,
287    ScaledIdentityOperator,
288    TFQMROptions,
289    TFQMRResult,
290};
291
292// Format conversions
293pub mod convert;
294
295// Construction utilities
296pub mod construct;
297pub mod construct_sym;
298
299// Combining arrays
300pub mod combine;
301pub use combine::{block_diag, bmat, hstack, kron, kronsum, tril, triu, vstack};
302
303// Index dtype handling utilities
304pub mod index_dtype;
305pub use index_dtype::{can_cast_safely, get_index_dtype, safely_cast_index_arrays};
306
307// Optimized operations for symmetric sparse formats
308pub mod sym_ops;
309
310// Tensor-based sparse operations (v0.2.0)
311pub mod tensor_sparse;
312
313// Enhanced BSR format with flat storage and Block LU (v0.3.0)
314pub mod bsr_enhanced;
315pub use bsr_enhanced::{block_lu, block_lu_solve, BlockLUResult, EnhancedBsr};
316
317// Enhanced DIA format with banded operations (v0.3.0)
318pub mod dia_enhanced;
319pub use dia_enhanced::{banded_solve, dia_tridiagonal_solve, tridiagonal_solve, EnhancedDia};
320
321// Compressed Sparse Fiber (CSF) tensor format (v0.3.0)
322pub mod csf_tensor;
323pub use csf_tensor::CsfTensor;
324
325// Sparse matrix utility functions (v0.3.0)
326pub mod sparse_functions;
327pub use sparse_functions::{
328    sparse_block_diag, sparse_diag_matrix, sparse_diags, sparse_eye, sparse_eye_rect,
329    sparse_hstack, sparse_kron, sparse_kronsum, sparse_random, sparse_vstack,
330};
331pub use sym_ops::{
332    sym_coo_matvec, sym_csr_matvec, sym_csr_quadratic_form, sym_csr_rank1_update, sym_csr_trace,
333};
334
335// Tensor operations (v0.2.0)
336pub use tensor_sparse::{khatri_rao_product, CPDecomposition, SparseTensor, TuckerDecomposition};
337
338// GPU-accelerated operations
339pub mod gpu;
340pub mod gpu_kernel_execution;
341pub mod gpu_ops;
342pub mod gpu_spmv_implementation;
343pub use gpu_kernel_execution::{
344    calculate_adaptive_workgroup_size, execute_spmv_kernel, execute_symmetric_spmv_kernel,
345    execute_triangular_solve_kernel, GpuKernelConfig, GpuMemoryManager as GpuKernelMemoryManager,
346    GpuPerformanceProfiler, MemoryStrategy,
347};
348pub use gpu_ops::{
349    gpu_sparse_matvec, gpu_sym_sparse_matvec, AdvancedGpuOps, GpuKernelScheduler, GpuMemoryManager,
350    GpuOptions, GpuProfiler, OptimizedGpuOps,
351};
352pub use gpu_spmv_implementation::GpuSpMV;
353
354// Memory-efficient algorithms and patterns
355pub mod memory_efficient;
356pub use memory_efficient::{
357    streaming_sparse_matvec, CacheAwareOps, MemoryPool, MemoryTracker, OutOfCoreProcessor,
358};
359
360// SIMD-accelerated operations
361pub mod simd_ops;
362pub use simd_ops::{
363    simd_csr_matvec, simd_sparse_elementwise, simd_sparse_linear_combination, simd_sparse_matmul,
364    simd_sparse_norm, simd_sparse_scale, simd_sparse_transpose, ElementwiseOp, SimdOptions,
365};
366
367// Parallel vector operations for iterative solvers
368pub mod parallel_vector_ops;
369pub use parallel_vector_ops::{
370    advanced_sparse_matvec_csr, parallel_axpy, parallel_dot, parallel_linear_combination,
371    parallel_norm2, parallel_sparse_matvec_csr, parallel_vector_add, parallel_vector_copy,
372    parallel_vector_scale, parallel_vector_sub, ParallelVectorOptions,
373};
374
375// Enhanced iterative solvers with preconditioners and sparse utilities (v0.3.0)
376pub mod iterative_solvers;
377pub use iterative_solvers::{
378    // Solvers
379    bicgstab as enhanced_bicgstab,
380    cg as enhanced_cg,
381    chebyshev,
382    // Sparse utility functions
383    estimate_spectral_radius,
384    gmres as enhanced_gmres,
385    sparse_diagonal,
386    sparse_norm,
387    sparse_trace,
388    // Preconditioners (Array1-based interface)
389    ILU0Preconditioner as EnhancedILU0Preconditioner,
390    // Configuration and result types
391    IterativeSolverConfig,
392    JacobiPreconditioner as EnhancedJacobiPreconditioner,
393    NormType,
394    Preconditioner,
395    SSORPreconditioner as EnhancedSSORPreconditioner,
396    SolverResult,
397};
398
399// LOBPCG eigensolver (v0.3.0)
400pub mod lobpcg;
401pub use lobpcg::{
402    lobpcg as lobpcg_eigensolver, lobpcg_generalised, EigenTarget, LobpcgConfig, LobpcgResult,
403};
404
405// Advanced Krylov subspace eigensolvers (v0.3.0)
406pub mod krylov;
407pub use krylov::{
408    iram as krylov_iram, thick_restart_lanczos, IramConfig, KrylovEigenResult,
409    ThickRestartLanczosConfig, WhichEigenvalues,
410};
411
412// Sparse matrix utilities (v0.3.0)
413pub mod sparse_utils;
414pub use sparse_utils::{
415    condest_1norm, permute_matrix, reverse_cuthill_mckee, sparse_add, sparse_extract_diagonal,
416    sparse_identity, sparse_kronecker, sparse_matrix_norm, sparse_matrix_trace, sparse_scale,
417    sparse_sub, sparse_transpose, spgemm, RcmResult, SparseNorm,
418};
419
420// Incomplete factorization preconditioners (v0.3.0)
421pub mod incomplete_factorizations;
422pub use incomplete_factorizations::{Ic0, Ilu0 as Ilu0Enhanced, IluK, Ilut, IlutConfig, Milu};
423
424// Sparse direct solvers (v0.3.0)
425pub mod direct_solver;
426pub use direct_solver::{
427    amd_ordering, elimination_tree, inverse_perm, nested_dissection_ordering,
428    sparse_cholesky_solve, sparse_lu_solve, symbolic_cholesky, SparseCholResult,
429    SparseCholeskySolver, SparseLuResult, SparseLuSolver, SparseSolver, SymbolicAnalysis,
430};
431
432// Sparse QR factorization (v0.3.0)
433pub mod sparse_qr;
434pub use sparse_qr::{
435    apply_q, apply_qt, extract_q_dense, numerical_rank, sparse_least_squares,
436    sparse_qr as sparse_qr_factorize, SparseLeastSquaresResult, SparseQrConfig, SparseQrResult,
437};
438
439// Unified sparse eigenvalue interface (v0.3.0)
440pub mod sparse_eigen;
441pub use sparse_eigen::{
442    cayley_transform_matvec, check_eigenpairs, compute_residuals, shift_invert_eig, sparse_eig,
443    sparse_eigs, sparse_eigsh, EigenMethod, EigenvalueTarget, SparseEigenConfig, SparseEigenResult,
444    SpectralTransform,
445};
446
447// Quantum-inspired sparse matrix operations (Advanced mode)
448pub mod quantum_inspired_sparse;
449pub use quantum_inspired_sparse::{
450    QuantumProcessorStats, QuantumSparseConfig, QuantumSparseProcessor, QuantumStrategy,
451};
452
453// Neural-adaptive sparse matrix operations (Advanced mode)
454pub mod neural_adaptive_sparse;
455pub use neural_adaptive_sparse::{
456    NeuralAdaptiveConfig, NeuralAdaptiveSparseProcessor, NeuralProcessorStats, OptimizationStrategy,
457};
458
459// Quantum-Neural hybrid optimization (Advanced mode)
460pub mod quantum_neural_hybrid;
461pub use quantum_neural_hybrid::{
462    HybridStrategy, QuantumNeuralConfig, QuantumNeuralHybridProcessor, QuantumNeuralHybridStats,
463};
464
465// Adaptive memory compression for advanced-large sparse matrices (Advanced mode)
466pub mod adaptive_memory_compression;
467pub use adaptive_memory_compression::{
468    AdaptiveCompressionConfig, AdaptiveMemoryCompressor, CompressedMatrix, CompressionAlgorithm,
469    MemoryStats,
470};
471
472// Real-time performance monitoring and adaptation (Advanced mode)
473pub mod realtime_performance_monitor;
474pub use realtime_performance_monitor::{
475    Alert, AlertSeverity, PerformanceMonitorConfig, PerformanceSample, ProcessorType,
476    RealTimePerformanceMonitor,
477};
478
479// Compressed sparse graph algorithms
480pub mod csgraph;
481pub use csgraph::{
482    all_pairs_shortest_path,
483    bellman_ford_single_source,
484    // Centrality measures (v0.2.0)
485    betweenness_centrality,
486    bfs_distances,
487    // Traversal algorithms
488    breadth_first_search,
489    closeness_centrality,
490    // Community detection (v0.2.0)
491    community_detection,
492    compute_laplacianmatrix,
493    connected_components,
494    degree_matrix,
495    depth_first_search,
496    dijkstra_single_source,
497    // Max flow algorithms (v0.2.0)
498    dinic,
499    edmonds_karp,
500    eigenvector_centrality,
501    floyd_warshall,
502    ford_fulkerson,
503    has_path,
504    is_connected,
505    is_laplacian,
506    is_spanning_tree,
507    // Minimum spanning trees
508    kruskal_mst,
509    label_propagation,
510    // Laplacian matrices
511    laplacian,
512    largest_component,
513    louvain_communities,
514    min_cut,
515    minimum_spanning_tree,
516    modularity,
517    num_edges,
518    num_vertices,
519    pagerank,
520    prim_mst,
521    reachable_vertices,
522    reconstruct_path,
523    // Graph algorithms
524    shortest_path,
525    // Shortest path algorithms
526    single_source_shortest_path,
527    spanning_tree_weight,
528    strongly_connected_components,
529    to_adjacency_list,
530    topological_sort,
531    traversegraph,
532    // Connected components
533    undirected_connected_components,
534    // Graph utilities
535    validate_graph,
536    weakly_connected_components,
537    LaplacianType,
538    MSTAlgorithm,
539    // Max flow types (v0.2.0)
540    MaxFlowResult,
541    // Enums and types
542    ShortestPathMethod,
543    TraversalOrder,
544};
545
546// Re-export warnings from scipy for compatibility
547pub struct SparseEfficiencyWarning;
548pub struct SparseWarning;
549
550/// Check if an object is a sparse array
551#[allow(dead_code)]
552pub fn is_sparse_array<T>(obj: &dyn SparseArray<T>) -> bool
553where
554    T: scirs2_core::SparseElement + std::ops::Div<Output = T> + PartialOrd + 'static,
555{
556    sparray::is_sparse(obj)
557}
558
559/// Check if an object is a symmetric sparse array
560#[allow(dead_code)]
561pub fn is_sym_sparse_array<T>(obj: &dyn SymSparseArray<T>) -> bool
562where
563    T: scirs2_core::SparseElement
564        + std::ops::Div<Output = T>
565        + scirs2_core::Float
566        + PartialOrd
567        + 'static,
568{
569    obj.is_symmetric()
570}
571
572/// Check if an object is a sparse matrix (legacy API)
573#[allow(dead_code)]
574pub fn is_sparse_matrix(obj: &dyn std::any::Any) -> bool {
575    obj.is::<CsrMatrix<f64>>()
576        || obj.is::<CscMatrix<f64>>()
577        || obj.is::<CooMatrix<f64>>()
578        || obj.is::<DokMatrix<f64>>()
579        || obj.is::<LilMatrix<f64>>()
580        || obj.is::<DiaMatrix<f64>>()
581        || obj.is::<BsrMatrix<f64>>()
582        || obj.is::<SymCsrMatrix<f64>>()
583        || obj.is::<SymCooMatrix<f64>>()
584        || obj.is::<CsrMatrix<f32>>()
585        || obj.is::<CscMatrix<f32>>()
586        || obj.is::<CooMatrix<f32>>()
587        || obj.is::<DokMatrix<f32>>()
588        || obj.is::<LilMatrix<f32>>()
589        || obj.is::<DiaMatrix<f32>>()
590        || obj.is::<BsrMatrix<f32>>()
591        || obj.is::<SymCsrMatrix<f32>>()
592        || obj.is::<SymCooMatrix<f32>>()
593}
594
595#[cfg(test)]
596mod tests {
597    use super::*;
598    use approx::assert_relative_eq;
599
600    #[test]
601    fn test_csr_array() {
602        let rows = vec![0, 0, 1, 2, 2];
603        let cols = vec![0, 2, 2, 0, 1];
604        let data = vec![1.0, 2.0, 3.0, 4.0, 5.0];
605        let shape = (3, 3);
606
607        let array =
608            CsrArray::from_triplets(&rows, &cols, &data, shape, false).expect("Operation failed");
609
610        assert_eq!(array.shape(), (3, 3));
611        assert_eq!(array.nnz(), 5);
612        assert!(is_sparse_array(&array));
613    }
614
615    #[test]
616    fn test_coo_array() {
617        let rows = vec![0, 0, 1, 2, 2];
618        let cols = vec![0, 2, 2, 0, 1];
619        let data = vec![1.0, 2.0, 3.0, 4.0, 5.0];
620        let shape = (3, 3);
621
622        let array =
623            CooArray::from_triplets(&rows, &cols, &data, shape, false).expect("Operation failed");
624
625        assert_eq!(array.shape(), (3, 3));
626        assert_eq!(array.nnz(), 5);
627        assert!(is_sparse_array(&array));
628    }
629
630    #[test]
631    fn test_dok_array() {
632        let rows = vec![0, 0, 1, 2, 2];
633        let cols = vec![0, 2, 2, 0, 1];
634        let data = vec![1.0, 2.0, 3.0, 4.0, 5.0];
635        let shape = (3, 3);
636
637        let array = DokArray::from_triplets(&rows, &cols, &data, shape).expect("Operation failed");
638
639        assert_eq!(array.shape(), (3, 3));
640        assert_eq!(array.nnz(), 5);
641        assert!(is_sparse_array(&array));
642
643        // Test setting and getting values
644        let mut array = DokArray::<f64>::new((2, 2));
645        array.set(0, 0, 1.0).expect("Operation failed");
646        array.set(1, 1, 2.0).expect("Operation failed");
647
648        assert_eq!(array.get(0, 0), 1.0);
649        assert_eq!(array.get(0, 1), 0.0);
650        assert_eq!(array.get(1, 1), 2.0);
651
652        // Test removing zeros
653        array.set(0, 0, 0.0).expect("Operation failed");
654        assert_eq!(array.nnz(), 1);
655    }
656
657    #[test]
658    fn test_lil_array() {
659        let rows = vec![0, 0, 1, 2, 2];
660        let cols = vec![0, 2, 2, 0, 1];
661        let data = vec![1.0, 2.0, 3.0, 4.0, 5.0];
662        let shape = (3, 3);
663
664        let array = LilArray::from_triplets(&rows, &cols, &data, shape).expect("Operation failed");
665
666        assert_eq!(array.shape(), (3, 3));
667        assert_eq!(array.nnz(), 5);
668        assert!(is_sparse_array(&array));
669
670        // Test setting and getting values
671        let mut array = LilArray::<f64>::new((2, 2));
672        array.set(0, 0, 1.0).expect("Operation failed");
673        array.set(1, 1, 2.0).expect("Operation failed");
674
675        assert_eq!(array.get(0, 0), 1.0);
676        assert_eq!(array.get(0, 1), 0.0);
677        assert_eq!(array.get(1, 1), 2.0);
678
679        // Test sorted indices
680        assert!(array.has_sorted_indices());
681
682        // Test removing zeros
683        array.set(0, 0, 0.0).expect("Operation failed");
684        assert_eq!(array.nnz(), 1);
685    }
686
687    #[test]
688    fn test_dia_array() {
689        use scirs2_core::ndarray::Array1;
690
691        // Create a 3x3 diagonal matrix with main diagonal + upper diagonal
692        let data = vec![
693            Array1::from_vec(vec![1.0, 2.0, 3.0]), // Main diagonal
694            Array1::from_vec(vec![4.0, 5.0, 0.0]), // Upper diagonal
695        ];
696        let offsets = vec![0, 1]; // Main diagonal and k=1
697        let shape = (3, 3);
698
699        let array = DiaArray::new(data, offsets, shape).expect("Operation failed");
700
701        assert_eq!(array.shape(), (3, 3));
702        assert_eq!(array.nnz(), 5); // 3 on main diagonal, 2 on upper diagonal
703        assert!(is_sparse_array(&array));
704
705        // Test values
706        assert_eq!(array.get(0, 0), 1.0);
707        assert_eq!(array.get(1, 1), 2.0);
708        assert_eq!(array.get(2, 2), 3.0);
709        assert_eq!(array.get(0, 1), 4.0);
710        assert_eq!(array.get(1, 2), 5.0);
711        assert_eq!(array.get(0, 2), 0.0);
712
713        // Test from_triplets
714        let rows = vec![0, 0, 1, 1, 2];
715        let cols = vec![0, 1, 1, 2, 2];
716        let data_vec = vec![1.0, 4.0, 2.0, 5.0, 3.0];
717
718        let array2 =
719            DiaArray::from_triplets(&rows, &cols, &data_vec, shape).expect("Operation failed");
720
721        // Should have same values
722        assert_eq!(array2.get(0, 0), 1.0);
723        assert_eq!(array2.get(1, 1), 2.0);
724        assert_eq!(array2.get(2, 2), 3.0);
725        assert_eq!(array2.get(0, 1), 4.0);
726        assert_eq!(array2.get(1, 2), 5.0);
727
728        // Test conversion to other formats
729        let csr = array.to_csr().expect("Operation failed");
730        assert_eq!(csr.nnz(), 5);
731        assert_eq!(csr.get(0, 0), 1.0);
732        assert_eq!(csr.get(0, 1), 4.0);
733    }
734
735    #[test]
736    fn test_format_conversions() {
737        let rows = vec![0, 0, 1, 2, 2];
738        let cols = vec![0, 2, 1, 0, 2];
739        let data = vec![1.0, 2.0, 3.0, 4.0, 5.0];
740        let shape = (3, 3);
741
742        // Create a COO array
743        let coo =
744            CooArray::from_triplets(&rows, &cols, &data, shape, false).expect("Operation failed");
745
746        // Convert to CSR
747        let csr = coo.to_csr().expect("Operation failed");
748
749        // Check values are preserved
750        let coo_dense = coo.to_array();
751        let csr_dense = csr.to_array();
752
753        for i in 0..shape.0 {
754            for j in 0..shape.1 {
755                assert_relative_eq!(coo_dense[[i, j]], csr_dense[[i, j]]);
756            }
757        }
758    }
759
760    #[test]
761    fn test_dot_product() {
762        let rows = vec![0, 0, 1, 2, 2];
763        let cols = vec![0, 2, 1, 0, 2];
764        let data = vec![1.0, 2.0, 3.0, 4.0, 5.0];
765        let shape = (3, 3);
766
767        // Create arrays in different formats
768        let coo =
769            CooArray::from_triplets(&rows, &cols, &data, shape, false).expect("Operation failed");
770        let csr =
771            CsrArray::from_triplets(&rows, &cols, &data, shape, false).expect("Operation failed");
772
773        // Compute dot product (matrix multiplication)
774        let coo_result = coo.dot(&coo).expect("Operation failed");
775        let csr_result = csr.dot(&csr).expect("Operation failed");
776
777        // Check results match
778        let coo_dense = coo_result.to_array();
779        let csr_dense = csr_result.to_array();
780
781        for i in 0..shape.0 {
782            for j in 0..shape.1 {
783                assert_relative_eq!(coo_dense[[i, j]], csr_dense[[i, j]], epsilon = 1e-10);
784            }
785        }
786    }
787
788    #[test]
789    fn test_sym_csr_array() {
790        // Create a symmetric matrix
791        let data = vec![2.0, 1.0, 2.0, 3.0, 0.0, 3.0, 1.0];
792        let indices = vec![0, 0, 1, 2, 0, 1, 2];
793        let indptr = vec![0, 1, 3, 7];
794
795        let sym_matrix =
796            SymCsrMatrix::new(data, indptr, indices, (3, 3)).expect("Operation failed");
797        let sym_array = SymCsrArray::new(sym_matrix);
798
799        assert_eq!(sym_array.shape(), (3, 3));
800        assert!(is_sym_sparse_array(&sym_array));
801
802        // Check values
803        assert_eq!(SparseArray::get(&sym_array, 0, 0), 2.0);
804        assert_eq!(SparseArray::get(&sym_array, 0, 1), 1.0);
805        assert_eq!(SparseArray::get(&sym_array, 1, 0), 1.0); // Symmetric element
806        assert_eq!(SparseArray::get(&sym_array, 1, 2), 3.0);
807        assert_eq!(SparseArray::get(&sym_array, 2, 1), 3.0); // Symmetric element
808
809        // Convert to standard CSR
810        let csr = SymSparseArray::to_csr(&sym_array).expect("Operation failed");
811        assert_eq!(csr.nnz(), 10); // Full matrix with symmetric elements
812    }
813
814    #[test]
815    fn test_sym_coo_array() {
816        // Create a symmetric matrix in COO format
817        let data = vec![2.0, 1.0, 2.0, 3.0, 1.0];
818        let rows = vec![0, 1, 1, 2, 2];
819        let cols = vec![0, 0, 1, 1, 2];
820
821        let sym_matrix = SymCooMatrix::new(data, rows, cols, (3, 3)).expect("Operation failed");
822        let sym_array = SymCooArray::new(sym_matrix);
823
824        assert_eq!(sym_array.shape(), (3, 3));
825        assert!(is_sym_sparse_array(&sym_array));
826
827        // Check values
828        assert_eq!(SparseArray::get(&sym_array, 0, 0), 2.0);
829        assert_eq!(SparseArray::get(&sym_array, 0, 1), 1.0);
830        assert_eq!(SparseArray::get(&sym_array, 1, 0), 1.0); // Symmetric element
831        assert_eq!(SparseArray::get(&sym_array, 1, 2), 3.0);
832        assert_eq!(SparseArray::get(&sym_array, 2, 1), 3.0); // Symmetric element
833
834        // Test from_triplets with enforce symmetry
835        // Input is intentionally asymmetric - will be fixed by enforce_symmetric=true
836        let rows2 = vec![0, 0, 1, 1, 2, 1, 0];
837        let cols2 = vec![0, 1, 1, 2, 2, 0, 2];
838        let data2 = vec![2.0, 1.5, 2.0, 3.5, 1.0, 0.5, 0.0];
839
840        let sym_array2 = SymCooArray::from_triplets(&rows2, &cols2, &data2, (3, 3), true)
841            .expect("Operation failed");
842
843        // Should average the asymmetric values
844        assert_eq!(SparseArray::get(&sym_array2, 0, 1), 1.0); // Average of 1.5 and 0.5
845        assert_eq!(SparseArray::get(&sym_array2, 1, 0), 1.0); // Symmetric element
846        assert_eq!(SparseArray::get(&sym_array2, 0, 2), 0.0); // Zero element
847    }
848
849    #[test]
850    fn test_construct_sym_utils() {
851        // Test creating an identity matrix
852        let eye = construct_sym::eye_sym_array::<f64>(3, "csr").expect("Operation failed");
853
854        assert_eq!(eye.shape(), (3, 3));
855        assert_eq!(SparseArray::get(&*eye, 0, 0), 1.0);
856        assert_eq!(SparseArray::get(&*eye, 1, 1), 1.0);
857        assert_eq!(SparseArray::get(&*eye, 2, 2), 1.0);
858        assert_eq!(SparseArray::get(&*eye, 0, 1), 0.0);
859
860        // Test creating a tridiagonal matrix - with coo format since csr had issues
861        let diag = vec![2.0, 2.0, 2.0];
862        let offdiag = vec![1.0, 1.0];
863
864        let tri =
865            construct_sym::tridiagonal_sym_array(&diag, &offdiag, "coo").expect("Operation failed");
866
867        assert_eq!(tri.shape(), (3, 3));
868        assert_eq!(SparseArray::get(&*tri, 0, 0), 2.0); // Main diagonal
869        assert_eq!(SparseArray::get(&*tri, 1, 1), 2.0);
870        assert_eq!(SparseArray::get(&*tri, 2, 2), 2.0);
871        assert_eq!(SparseArray::get(&*tri, 0, 1), 1.0); // Off-diagonal
872        assert_eq!(SparseArray::get(&*tri, 1, 0), 1.0); // Symmetric element
873        assert_eq!(SparseArray::get(&*tri, 1, 2), 1.0);
874        assert_eq!(SparseArray::get(&*tri, 0, 2), 0.0); // Zero element
875
876        // Test creating a banded matrix
877        let diagonals = vec![
878            vec![2.0, 2.0, 2.0, 2.0, 2.0], // Main diagonal
879            vec![1.0, 1.0, 1.0, 1.0],      // First off-diagonal
880            vec![0.5, 0.5, 0.5],           // Second off-diagonal
881        ];
882
883        let band = construct_sym::banded_sym_array(&diagonals, 5, "csr").expect("Operation failed");
884
885        assert_eq!(band.shape(), (5, 5));
886        assert_eq!(SparseArray::get(&*band, 0, 0), 2.0);
887        assert_eq!(SparseArray::get(&*band, 0, 1), 1.0);
888        assert_eq!(SparseArray::get(&*band, 0, 2), 0.5);
889        assert_eq!(SparseArray::get(&*band, 2, 0), 0.5); // Symmetric element
890    }
891
892    #[test]
893    fn test_sym_conversions() {
894        // Create a symmetric matrix
895        // Lower triangular part only
896        let data = vec![2.0, 1.0, 2.0, 3.0, 1.0];
897        let rows = vec![0, 1, 1, 2, 2];
898        let cols = vec![0, 0, 1, 1, 2];
899
900        let sym_coo = SymCooArray::from_triplets(&rows, &cols, &data, (3, 3), true)
901            .expect("Operation failed");
902
903        // Convert to symmetric CSR
904        let sym_csr = sym_coo.to_sym_csr().expect("Operation failed");
905
906        // Check values are preserved
907        for i in 0..3 {
908            for j in 0..3 {
909                assert_eq!(
910                    SparseArray::get(&sym_coo, i, j),
911                    SparseArray::get(&sym_csr, i, j)
912                );
913            }
914        }
915
916        // Convert to standard formats
917        let csr = SymSparseArray::to_csr(&sym_coo).expect("Operation failed");
918        let coo = SymSparseArray::to_coo(&sym_csr).expect("Operation failed");
919
920        // Check full symmetric matrix in standard formats
921        assert_eq!(csr.nnz(), 7); // Accounts for symmetric pairs
922        assert_eq!(coo.nnz(), 7);
923
924        for i in 0..3 {
925            for j in 0..3 {
926                assert_eq!(SparseArray::get(&csr, i, j), SparseArray::get(&coo, i, j));
927                assert_eq!(
928                    SparseArray::get(&csr, i, j),
929                    SparseArray::get(&sym_csr, i, j)
930                );
931            }
932        }
933    }
934}