scirs2_linalg/matrix_functions/
analysis.rs

1//! Matrix analysis functions including decompositions and matrix norms
2
3use scirs2_core::ndarray::{Array2, ArrayView2};
4use scirs2_core::numeric::{Float, NumAssign, One};
5use std::iter::Sum;
6
7use crate::eigen::eig;
8use crate::error::{LinalgError, LinalgResult};
9use crate::solve::solve_multiple;
10use crate::validation::validate_decomposition;
11
12/// Compute the spectral radius of a matrix.
13///
14/// The spectral radius is the largest absolute value of the eigenvalues.
15///
16/// # Arguments
17///
18/// * `a` - Input square matrix
19/// * `workers` - Number of worker threads (None = use default)
20///
21/// # Returns
22///
23/// * Spectral radius (largest eigenvalue magnitude)
24///
25/// # Examples
26///
27/// ```
28/// use scirs2_core::ndarray::array;
29/// use scirs2_linalg::matrix_functions::spectral_radius;
30///
31/// let a = array![[2.0_f64, 0.0], [0.0, 3.0]];
32/// let rho = spectral_radius(&a.view(), None).unwrap();
33/// // rho should be 3.0
34/// ```
35#[allow(dead_code)]
36pub fn spectral_radius<F>(a: &ArrayView2<F>, workers: Option<usize>) -> LinalgResult<F>
37where
38    F: Float + NumAssign + Sum + One + Send + Sync + scirs2_core::ndarray::ScalarOperand + 'static,
39{
40    use crate::parallel;
41
42    // Configure workers for parallel operations
43    parallel::configure_workers(workers);
44
45    validate_decomposition(a, "Spectral radius computation", true)?;
46
47    // Compute eigenvalues
48    let (eigenvals, _) = eig(a, None)?;
49
50    // Find the maximum absolute value
51    let mut max_abs = F::zero();
52    for &val in eigenvals.iter() {
53        let abs_val = (val.re * val.re + val.im * val.im).sqrt();
54        if abs_val > max_abs {
55            max_abs = abs_val;
56        }
57    }
58
59    Ok(max_abs)
60}
61
62/// Compute the spectral condition number of a matrix.
63///
64/// The spectral condition number is the ratio of the largest to smallest
65/// singular values (or eigenvalues for symmetric matrices).
66///
67/// # Arguments
68///
69/// * `a` - Input square matrix
70/// * `workers` - Number of worker threads (None = use default)
71///
72/// # Returns
73///
74/// * Spectral condition number
75///
76/// # Examples
77///
78/// ```
79/// use scirs2_core::ndarray::array;
80/// use scirs2_linalg::matrix_functions::spectral_condition_number;
81///
82/// let a = array![[2.0_f64, 0.0], [0.0, 1.0]];
83/// let cond = spectral_condition_number(&a.view(), None).unwrap();
84/// // cond should be 2.0
85/// ```
86#[allow(dead_code)]
87pub fn spectral_condition_number<F>(a: &ArrayView2<F>, workers: Option<usize>) -> LinalgResult<F>
88where
89    F: Float + NumAssign + Sum + One + Send + Sync + scirs2_core::ndarray::ScalarOperand + 'static,
90{
91    use crate::parallel;
92
93    // Configure workers for parallel operations
94    parallel::configure_workers(workers);
95
96    validate_decomposition(a, "Spectral condition number computation", true)?;
97
98    // Compute eigenvalues
99    let (eigenvals, _) = eig(a, None)?;
100
101    // Find the maximum and minimum absolute values
102    let mut max_abs = F::zero();
103    let mut min_abs = F::infinity();
104
105    for &val in eigenvals.iter() {
106        let abs_val = (val.re * val.re + val.im * val.im).sqrt();
107        if abs_val > max_abs {
108            max_abs = abs_val;
109        }
110        if abs_val < min_abs && abs_val > F::epsilon() {
111            min_abs = abs_val;
112        }
113    }
114
115    // Check for singular matrix
116    if min_abs < F::epsilon() {
117        return Ok(F::infinity());
118    }
119
120    Ok(max_abs / min_abs)
121}
122
123/// Compute the polar decomposition of a matrix.
124///
125/// The polar decomposition factorizes a matrix A as A = UP where U is unitary
126/// and P is positive semidefinite.
127///
128/// # Arguments
129///
130/// * `a` - Input matrix
131/// * `side` - Which factor to return ("left" for A = UP, "right" for A = PU)
132///
133/// # Returns
134///
135/// * (U, P) - Unitary and positive semidefinite factors
136///
137/// # Examples
138///
139/// ```
140/// use scirs2_core::ndarray::array;
141/// use scirs2_linalg::matrix_functions::polar_decomposition;
142///
143/// let a = array![[2.0_f64, 1.0], [0.0, 1.0]];
144/// let (u, p) = polar_decomposition(&a.view(), "right").unwrap();
145/// ```
146#[allow(dead_code)]
147pub fn polar_decomposition<F>(a: &ArrayView2<F>, side: &str) -> LinalgResult<(Array2<F>, Array2<F>)>
148where
149    F: Float + NumAssign + Sum + One + Send + Sync + scirs2_core::ndarray::ScalarOperand + 'static,
150{
151    use super::exponential::sqrtm;
152
153    validate_decomposition(a, "Polar decomposition", false)?;
154
155    let (m, n) = a.dim();
156
157    match side {
158        "right" => {
159            // Right polar decomposition: A = UP
160            // P = sqrt(A^H * A), U = A * P^{-1}
161
162            // Compute A^H * A
163            let mut aha = Array2::<F>::zeros((n, n));
164            for i in 0..n {
165                for j in 0..n {
166                    for k in 0..m {
167                        aha[[i, j]] += a[[k, i]] * a[[k, j]]; // A^H * A
168                    }
169                }
170            }
171
172            // Compute P = sqrt(A^H * A)
173            let p = sqrtm(&aha.view(), 50, F::from(1e-12).unwrap())?;
174
175            // Compute U = A * P^{-1}
176            let p_inv = solve_multiple(&p.view(), &Array2::eye(n).view(), None)?;
177
178            let mut u = Array2::<F>::zeros((m, n));
179            for i in 0..m {
180                for j in 0..n {
181                    for k in 0..n {
182                        u[[i, j]] += a[[i, k]] * p_inv[[k, j]];
183                    }
184                }
185            }
186
187            Ok((u, p))
188        }
189        "left" => {
190            // Left polar decomposition: A = PU
191            // P = sqrt(A * A^H), U = P^{-1} * A
192
193            // Compute A * A^H
194            let mut aah = Array2::<F>::zeros((m, m));
195            for i in 0..m {
196                for j in 0..m {
197                    for k in 0..n {
198                        aah[[i, j]] += a[[i, k]] * a[[j, k]]; // A * A^H
199                    }
200                }
201            }
202
203            // Compute P = sqrt(A * A^H)
204            let p = sqrtm(&aah.view(), 50, F::from(1e-12).unwrap())?;
205
206            // Compute U = P^{-1} * A
207            let p_inv = solve_multiple(&p.view(), &Array2::eye(m).view(), None)?;
208
209            let mut u = Array2::<F>::zeros((m, n));
210            for i in 0..m {
211                for j in 0..n {
212                    for k in 0..m {
213                        u[[i, j]] += p_inv[[i, k]] * a[[k, j]];
214                    }
215                }
216            }
217
218            Ok((p, u))
219        }
220        _ => Err(LinalgError::InvalidInputError(format!(
221            "Invalid side '{}'. Must be 'left' or 'right'",
222            side
223        ))),
224    }
225}
226
227/// Compute the geometric mean of symmetric positive definite matrices.
228///
229/// For two SPD matrices A and B, the geometric mean is:
230/// A #_t B = A^{1/2} * (A^{-1/2} * B * A^{-1/2})^t * A^{1/2}
231///
232/// For t = 1/2, this gives the standard geometric mean.
233///
234/// # Arguments
235///
236/// * `a` - First SPD matrix
237/// * `b` - Second SPD matrix
238/// * `t` - Parameter (0.5 for standard geometric mean)
239///
240/// # Returns
241///
242/// * Geometric mean of A and B
243///
244/// # Examples
245///
246/// ```
247/// use scirs2_core::ndarray::array;
248/// use scirs2_linalg::matrix_functions::geometric_mean_spd;
249///
250/// let a = array![[4.0_f64, 0.0], [0.0, 1.0]];
251/// let b = array![[1.0_f64, 0.0], [0.0, 4.0]];
252/// let mean = geometric_mean_spd(&a.view(), &b.view(), 0.5).unwrap();
253/// ```
254#[allow(dead_code)]
255pub fn geometric_mean_spd<F>(a: &ArrayView2<F>, b: &ArrayView2<F>, t: F) -> LinalgResult<Array2<F>>
256where
257    F: Float + NumAssign + Sum + One + Send + Sync + scirs2_core::ndarray::ScalarOperand + 'static,
258{
259    use super::exponential::sqrtm;
260    use super::fractional::spdmatrix_function;
261
262    validate_decomposition(a, "Geometric mean computation (matrix A)", true)?;
263    validate_decomposition(b, "Geometric mean computation (matrix B)", true)?;
264
265    if a.dim() != b.dim() {
266        return Err(LinalgError::ShapeError(
267            "Matrices must have the same dimensions".to_string(),
268        ));
269    }
270
271    let n = a.nrows();
272
273    // Check if matrices are symmetric
274    for i in 0..n {
275        for j in 0..n {
276            if (a[[i, j]] - a[[j, i]]).abs() > F::epsilon() {
277                return Err(LinalgError::InvalidInputError(
278                    "Matrix A must be symmetric".to_string(),
279                ));
280            }
281            if (b[[i, j]] - b[[j, i]]).abs() > F::epsilon() {
282                return Err(LinalgError::InvalidInputError(
283                    "Matrix B must be symmetric".to_string(),
284                ));
285            }
286        }
287    }
288
289    // Compute A^{1/2}
290    let a_half = spdmatrix_function(a, |x| x.sqrt(), true)?;
291
292    // Compute A^{-1/2}
293    let a_neg_half = spdmatrix_function(a, |x| F::one() / x.sqrt(), true)?;
294
295    // Compute A^{-1/2} * B * A^{-1/2}
296    let mut temp1 = Array2::<F>::zeros((n, n));
297    for i in 0..n {
298        for j in 0..n {
299            for k in 0..n {
300                temp1[[i, j]] += a_neg_half[[i, k]] * b[[k, j]];
301            }
302        }
303    }
304
305    let mut similarity = Array2::<F>::zeros((n, n));
306    for i in 0..n {
307        for j in 0..n {
308            for k in 0..n {
309                similarity[[i, j]] += temp1[[i, k]] * a_neg_half[[k, j]];
310            }
311        }
312    }
313
314    // Compute (A^{-1/2} * B * A^{-1/2})^t
315    let powered_similarity = spdmatrix_function(&similarity.view(), |x| x.powf(t), true)?;
316
317    // Compute A^{1/2} * (A^{-1/2} * B * A^{-1/2})^t * A^{1/2}
318    let mut temp2 = Array2::<F>::zeros((n, n));
319    for i in 0..n {
320        for j in 0..n {
321            for k in 0..n {
322                temp2[[i, j]] += a_half[[i, k]] * powered_similarity[[k, j]];
323            }
324        }
325    }
326
327    let mut result = Array2::<F>::zeros((n, n));
328    for i in 0..n {
329        for j in 0..n {
330            for k in 0..n {
331                result[[i, j]] += temp2[[i, k]] * a_half[[k, j]];
332            }
333        }
334    }
335
336    Ok(result)
337}
338
339/// Compute Tikhonov regularization of a matrix.
340///
341/// Tikhonov regularization adds a ridge term to improve conditioning:
342/// A_reg = A + λI
343///
344/// # Arguments
345///
346/// * `a` - Input matrix
347/// * `lambda` - Regularization parameter
348/// * `identity_like` - Whether to add λI (true) or λ times identity-like term
349///
350/// # Returns
351///
352/// * Regularized matrix
353///
354/// # Examples
355///
356/// ```
357/// use scirs2_core::ndarray::array;
358/// use scirs2_linalg::matrix_functions::tikhonov_regularization;
359///
360/// let a = array![[1.0_f64, 0.5], [0.5, 1.0]];
361/// let reg_a = tikhonov_regularization(&a.view(), 0.1, true).unwrap();
362/// ```
363#[allow(dead_code)]
364pub fn tikhonov_regularization<F>(
365    a: &ArrayView2<F>,
366    lambda: F,
367    identity_like: bool,
368) -> LinalgResult<Array2<F>>
369where
370    F: Float + NumAssign + Sum + One + Send + Sync + scirs2_core::ndarray::ScalarOperand + 'static,
371{
372    let (m, n) = a.dim();
373
374    if identity_like && m != n {
375        return Err(LinalgError::ShapeError(
376            "Matrix must be square for identity-like regularization".to_string(),
377        ));
378    }
379
380    let mut result = a.to_owned();
381
382    if identity_like {
383        // Add λI to the diagonal
384        for i in 0..n {
385            result[[i, i]] += lambda;
386        }
387    } else {
388        // Add λ to all elements (different regularization scheme)
389        for i in 0..m {
390            for j in 0..n {
391                result[[i, j]] += lambda;
392            }
393        }
394    }
395
396    Ok(result)
397}
398
399/// Compute the nuclear norm (trace norm) of a matrix.
400///
401/// The nuclear norm is the sum of singular values, which for symmetric matrices
402/// is the sum of absolute eigenvalues.
403///
404/// # Arguments
405///
406/// * `a` - Input matrix
407/// * `workers` - Number of worker threads (None = use default)
408///
409/// # Returns
410///
411/// * Nuclear norm of the matrix
412///
413/// # Examples
414///
415/// ```
416/// use scirs2_core::ndarray::array;
417/// use scirs2_linalg::matrix_functions::nuclear_norm;
418///
419/// let a = array![[2.0_f64, 0.0], [0.0, 3.0]];
420/// let norm = nuclear_norm(&a.view(), None).unwrap();
421/// // norm should be 5.0
422/// ```
423#[allow(dead_code)]
424pub fn nuclear_norm<F>(a: &ArrayView2<F>, workers: Option<usize>) -> LinalgResult<F>
425where
426    F: Float + NumAssign + Sum + One + Send + Sync + scirs2_core::ndarray::ScalarOperand + 'static,
427{
428    use crate::parallel;
429
430    // Configure workers for parallel operations
431    parallel::configure_workers(workers);
432
433    validate_decomposition(a, "Nuclear norm computation", false)?;
434
435    let (m, n) = a.dim();
436
437    if m == n {
438        // For square matrices, use eigendecomposition
439        // This gives us the singular values for symmetric matrices
440        let (eigenvals, _) = eig(a, None)?;
441
442        let mut sum = F::zero();
443        for &val in eigenvals.iter() {
444            sum += (val.re * val.re + val.im * val.im).sqrt();
445        }
446
447        Ok(sum)
448    } else {
449        // For non-square matrices, we'd need SVD
450        // For now, return an error
451        Err(LinalgError::ImplementationError(
452            "Nuclear norm for non-square matrices requires SVD (not yet implemented)".to_string(),
453        ))
454    }
455}