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
//! Matrix exponential computation for sparse matrices
//!
//! This module implements the matrix exponential using the scaling and squaring
//! method with Padé approximation.

use crate::csr::CsrMatrix;
use crate::error::{SparseError, SparseResult};
use scirs2_core::numeric::{Float, NumAssign, One, SparseElement, Zero};
use std::iter::Sum;

/// Compute the matrix exponential using scaling and squaring with Padé approximation
///
/// This function computes exp(A) for a sparse matrix A using the scaling and
/// squaring method combined with Padé approximation.
///
/// # Arguments
///
/// * `a` - The sparse matrix A (must be square)
///
/// # Returns
///
/// The matrix exponential exp(A) as a sparse matrix
///
/// # Implementation Details
///
/// Uses 13th order Padé approximation for high accuracy (machine precision).
/// The algorithm automatically selects the appropriate scaling factor based
/// on the matrix norm to ensure numerical stability.
#[allow(dead_code)]
pub fn expm<F>(a: &CsrMatrix<F>) -> SparseResult<CsrMatrix<F>>
where
    F: Float + NumAssign + Sum + SparseElement + 'static + std::fmt::Debug,
{
    let (rows, cols) = a.shape();
    if rows != cols {
        return Err(SparseError::ValueError(
            "Matrix must be square for expm".to_string(),
        ));
    }

    // Compute the matrix infinity norm
    let a_norm = matrix_inf_norm(a)?;

    // Constants for order 13 Padé approximation
    let theta_13 = F::from(5.371920351148152).expect("Failed to convert constant to float");

    // If the norm is small enough, use direct Padé approximation
    if a_norm <= theta_13 {
        return pade_approximation(a, 13);
    }

    // Otherwise, use scaling and squaring
    // Find s such that ||A/2^s|| <= theta_13
    let mut s = 0;
    let mut scaled_norm = a_norm;
    let two = F::from(2.0).expect("Failed to convert constant to float");

    while scaled_norm > theta_13 {
        s += 1;
        scaled_norm /= two;
    }

    // Compute A/2^s
    let scale_factor = two.powi(s);
    let scaled_a = scale_matrix(a, F::sparse_one() / scale_factor)?;

    // Compute exp(A/2^s) using Padé approximation
    let mut exp_scaled = pade_approximation(&scaled_a, 13)?;

    // Square the result s times to get exp(A)
    for _ in 0..s {
        exp_scaled = exp_scaled.matmul(&exp_scaled)?;
    }

    Ok(exp_scaled)
}

/// Compute the Padé approximation of exp(A)
///
/// Uses the diagonal Padé approximant of order (p,p)
#[allow(dead_code)]
fn pade_approximation<F>(a: &CsrMatrix<F>, p: usize) -> SparseResult<CsrMatrix<F>>
where
    F: Float + NumAssign + Sum + SparseElement + 'static + std::fmt::Debug,
{
    let n = a.shape().0;

    // Compute powers of A
    let mut a_powers = vec![sparse_identity(n)?]; // A^0 = I
    a_powers.push(a.clone()); // A^1 = A

    // Compute A^2, A^3, ..., A^p
    for i in 2..=p {
        let prev = &a_powers[i - 1];
        let power = prev.matmul(a)?;
        a_powers.push(power);
    }

    // Compute Padé coefficients
    let pade_coeffs = match p {
        6 => vec![
            F::from(1.0).expect("Failed to convert constant to float"),
            F::from(1.0 / 2.0).expect("Failed to convert to float"),
            F::from(3.0 / 26.0).expect("Failed to convert to float"),
            F::from(1.0 / 312.0).expect("Failed to convert to float"),
            F::from(1.0 / 10608.0).expect("Failed to convert to float"),
            F::from(1.0 / 358800.0).expect("Failed to convert to float"),
            F::from(1.0 / 17297280.0).expect("Failed to convert to float"),
        ],
        13 => {
            // Compute coefficients for Padé (13,13) approximant
            // c_k = (2p-k)! p! / ((2p)! k! (p-k)!) for k = 0, 1, ..., p
            let two_p = 26i64;
            let p = 13i64;
            let mut coeffs = Vec::with_capacity(14);

            for k in 0..=p {
                let mut num = 1.0;
                let mut den = 1.0;

                // (2p-k)! / (2p)!
                for i in (two_p - k + 1)..=two_p {
                    den *= i as f64;
                }

                // p! / (p-k)!
                for i in (p - k + 1)..=p {
                    num *= i as f64;
                }

                // 1 / k!
                let mut k_fact = 1.0;
                for i in 1..=k {
                    k_fact *= i as f64;
                }

                coeffs.push(F::from(num / (den * k_fact)).expect("Operation failed"));
            }

            coeffs
        }
        _ => {
            // General formula for Padé coefficients
            let mut coeffs = vec![F::sparse_zero(); p + 1];
            let mut factorial: F = F::sparse_one();
            for (i, coeff) in coeffs.iter_mut().enumerate().take(p + 1) {
                if i > 0 {
                    factorial *= F::from(i).expect("Failed to convert to float");
                }
                let numerator = factorial;
                let mut denominator = F::sparse_one();
                for j in 1..=i {
                    denominator *= F::from(p + 1 - j).expect("Failed to convert to float");
                }
                for j in 1..=(p - i) {
                    denominator *= F::from(j).expect("Failed to convert to float");
                }
                *coeff = numerator / denominator;
            }
            coeffs
        }
    };

    // Compute U and V for the Padé approximant
    let mut u = sparse_zero(n)?;
    let mut v = sparse_zero(n)?;

    // U = sum of odd powers, V = sum of even powers
    for (i, coeff) in pade_coeffs.iter().enumerate() {
        let scaled_matrix = scale_matrix(&a_powers[i], *coeff)?;
        if i % 2 == 0 {
            v = sparse_add(&v, &scaled_matrix)?;
        } else {
            u = sparse_add(&u, &scaled_matrix)?;
        }
    }

    // Compute (V - U)^(-1) * (V + U)
    let neg_u = scale_matrix(
        &u,
        F::from(-1.0).expect("Failed to convert constant to float"),
    )?;
    let v_minus_u = sparse_add(&v, &neg_u)?;
    let v_plus_u = sparse_add(&v, &u)?;

    // Solve (V - U) * X = (V + U) for X
    sparse_solve(&v_minus_u, &v_plus_u)
}

/// Compute the infinity norm of a sparse matrix
#[allow(dead_code)]
fn matrix_inf_norm<F>(a: &CsrMatrix<F>) -> SparseResult<F>
where
    F: Float + NumAssign + Sum + SparseElement + std::fmt::Debug,
{
    let mut max_row_sum = F::sparse_zero();

    // For CSR format, efficiently compute row sums
    for row in 0..a.rows() {
        let start = a.indptr[row];
        let end = a.indptr[row + 1];
        let row_sum: F = a.data[start..end].iter().map(|x| x.abs()).sum();

        if row_sum > max_row_sum {
            max_row_sum = row_sum;
        }
    }

    Ok(max_row_sum)
}

/// Scale a sparse matrix by a scalar
#[allow(dead_code)]
fn scale_matrix<F>(a: &CsrMatrix<F>, scale: F) -> SparseResult<CsrMatrix<F>>
where
    F: Float + NumAssign + SparseElement,
{
    let mut data = a.data.clone();
    for val in data.iter_mut() {
        *val *= scale;
    }
    CsrMatrix::from_raw_csr(data, a.indptr.clone(), a.indices.clone(), a.shape())
}

/// Create a sparse identity matrix
#[allow(dead_code)]
fn sparse_identity<F>(n: usize) -> SparseResult<CsrMatrix<F>>
where
    F: Float + Zero + One + SparseElement,
{
    let mut rows = Vec::with_capacity(n);
    let mut cols = Vec::with_capacity(n);
    let mut values = Vec::with_capacity(n);

    for i in 0..n {
        rows.push(i);
        cols.push(i);
        values.push(F::sparse_one());
    }

    CsrMatrix::new(values, rows, cols, (n, n))
}

/// Create a sparse zero matrix
#[allow(dead_code)]
fn sparse_zero<F>(n: usize) -> SparseResult<CsrMatrix<F>>
where
    F: Float + Zero + SparseElement,
{
    Ok(CsrMatrix::empty((n, n)))
}

/// Add two sparse matrices
#[allow(dead_code)]
fn sparse_add<F>(a: &CsrMatrix<F>, b: &CsrMatrix<F>) -> SparseResult<CsrMatrix<F>>
where
    F: Float + NumAssign + SparseElement,
{
    if a.shape() != b.shape() {
        return Err(SparseError::ShapeMismatch {
            expected: a.shape(),
            found: b.shape(),
        });
    }

    let mut rows = Vec::new();
    let mut cols = Vec::new();
    let mut values = Vec::new();

    for i in 0..a.rows() {
        for j in 0..a.cols() {
            let val = a.get(i, j) + b.get(i, j);
            if val.abs() > F::epsilon() {
                rows.push(i);
                cols.push(j);
                values.push(val);
            }
        }
    }

    CsrMatrix::new(values, rows, cols, a.shape())
}

/// Solve a linear system A * X = B for sparse matrices
///
/// Uses BiCGSTAB iterative solver column-by-column to solve AX = B.
#[allow(dead_code)]
fn sparse_solve<F>(a: &CsrMatrix<F>, b: &CsrMatrix<F>) -> SparseResult<CsrMatrix<F>>
where
    F: Float + NumAssign + Sum + SparseElement + 'static + std::fmt::Debug,
{
    use crate::linalg::interface::MatrixLinearOperator;
    use crate::linalg::iterative::bicgstab;
    use crate::linalg::iterative::BiCGSTABOptions;

    let n = a.rows();
    let mut result_rows = Vec::new();
    let mut result_cols = Vec::new();
    let mut result_values = Vec::new();

    // Solve column by column
    for col in 0..b.cols() {
        // Extract the column from B
        let b_col = (0..n).map(|row| b.get(row, col)).collect::<Vec<_>>();

        // Create a linear operator for the matrix
        let op = MatrixLinearOperator::new(a.clone());

        // Create options for BiCGSTAB
        let options = BiCGSTABOptions {
            rtol: F::from(1e-10).expect("Failed to convert constant to float"),
            atol: F::from(1e-12).expect("Failed to convert constant to float"),
            max_iter: 1000,
            x0: None,
            left_preconditioner: None,
            right_preconditioner: None,
        };

        // Use BiCGSTAB to solve A * x = b_col
        let result = bicgstab(&op, &b_col, options)?;

        // Check convergence
        if !result.converged {
            return Err(SparseError::IterativeSolverFailure(format!(
                "BiCGSTAB failed to converge in {} iterations",
                result.iterations
            )));
        }

        // Add non-zero entries to result
        for (row, &val) in result.x.iter().enumerate() {
            if val.abs() > F::epsilon() {
                result_rows.push(row);
                result_cols.push(col);
                result_values.push(val);
            }
        }
    }

    CsrMatrix::new(result_values, result_rows, result_cols, (n, b.cols()))
}

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

    #[test]
    fn test_expm_identity() {
        // exp(0) = I
        let n = 3;
        let zero_matrix = sparse_zero::<f64>(n).expect("Operation failed");
        let exp_zero = expm(&zero_matrix).expect("Operation failed");

        // Check that exp(0) is identity
        for i in 0..n {
            for j in 0..n {
                let expected = if i == j { 1.0 } else { 0.0 };
                let actual = exp_zero.get(i, j);
                assert_relative_eq!(actual, expected, epsilon = 1e-10);
            }
        }
    }

    #[test]
    fn test_expm_diagonal() {
        // For diagonal matrix D, exp(D) is diagonal with exp(d_ii) on diagonal
        let n = 3;
        let diag_values = [0.5, 1.0, 2.0];
        let mut rows = Vec::new();
        let mut cols = Vec::new();
        let mut values = Vec::new();

        for (i, &val) in diag_values.iter().enumerate() {
            rows.push(i);
            cols.push(i);
            values.push(val);
        }

        let diag_matrix = CsrMatrix::new(values, rows, cols, (n, n)).expect("Operation failed");
        let exp_diag = expm(&diag_matrix).expect("Operation failed");

        // Check diagonal values with high precision
        for (i, &val) in diag_values.iter().enumerate() {
            let expected = val.exp();
            let actual = exp_diag.get(i, i);
            assert_relative_eq!(actual, expected, epsilon = 1e-10);
        }

        // Check off-diagonal values are zero
        for i in 0..n {
            for j in 0..n {
                if i != j {
                    let actual = exp_diag.get(i, j);
                    assert_relative_eq!(actual, 0.0, epsilon = 1e-10);
                }
            }
        }
    }

    #[test]
    fn test_expm_small_matrix() {
        // Test on a small 2x2 matrix with known exponential
        // A = [[0, 1], [0, 0]]
        // exp(A) = [[1, 1], [0, 1]]
        let rows = vec![0, 1];
        let cols = vec![1, 0];
        let values = vec![1.0, 0.0];

        let a = CsrMatrix::new(values, rows, cols, (2, 2)).expect("Operation failed");
        let exp_a = expm(&a).expect("Operation failed");

        // Check the result
        assert_relative_eq!(exp_a.get(0, 0), 1.0, epsilon = 1e-10);
        assert_relative_eq!(exp_a.get(0, 1), 1.0, epsilon = 1e-10);
        assert_relative_eq!(exp_a.get(1, 0), 0.0, epsilon = 1e-10);
        assert_relative_eq!(exp_a.get(1, 1), 1.0, epsilon = 1e-10);
    }
}