use scirs2_core::ndarray::Array2;
use scirs2_core::random::{Rng, SeedableRng, StdRng};
use crate::error::{GraphError, Result};
pub fn graph_laplacian(adj: &Array2<f64>) -> Array2<f64> {
let n = adj.nrows();
let mut lap = Array2::zeros((n, n));
for i in 0..n {
let deg: f64 = adj.row(i).sum();
lap[[i, i]] = deg;
for j in 0..n {
if i != j {
lap[[i, j]] = -adj[[i, j]];
}
}
}
lap
}
pub fn normalized_laplacian(adj: &Array2<f64>) -> Result<Array2<f64>> {
let n = adj.nrows();
if n == 0 {
return Err(GraphError::InvalidGraph("empty adjacency matrix".into()));
}
if adj.ncols() != n {
return Err(GraphError::InvalidGraph("adjacency matrix must be square".into()));
}
let d_inv_sqrt: Vec<f64> = (0..n)
.map(|i| {
let deg = adj.row(i).sum();
if deg > 0.0 {
1.0 / deg.sqrt()
} else {
0.0
}
})
.collect();
let l = graph_laplacian(adj);
let mut l_norm = Array2::zeros((n, n));
for i in 0..n {
for j in 0..n {
l_norm[[i, j]] = d_inv_sqrt[i] * l[[i, j]] * d_inv_sqrt[j];
}
}
Ok(l_norm)
}
fn matvec(m: &Array2<f64>, v: &[f64]) -> Vec<f64> {
let n = m.nrows();
let mut out = vec![0.0f64; n];
for i in 0..n {
for j in 0..n {
out[i] += m[[i, j]] * v[j];
}
}
out
}
fn vec_norm(v: &[f64]) -> f64 {
v.iter().map(|x| x * x).sum::<f64>().sqrt()
}
fn normalize_inplace(v: &mut Vec<f64>) -> f64 {
let n = vec_norm(v);
if n > 1e-300 {
v.iter_mut().for_each(|x| *x /= n);
}
n
}
fn deflate(v: &mut Vec<f64>, basis: &[Vec<f64>]) {
for b in basis {
let dot: f64 = v.iter().zip(b.iter()).map(|(a, c)| a * c).sum();
for (vi, bi) in v.iter_mut().zip(b.iter()) {
*vi -= dot * bi;
}
}
}
fn power_iteration(m: &Array2<f64>, max_iter: usize, tol: f64, seed: u64) -> (f64, Vec<f64>) {
let n = m.nrows();
let mut rng = StdRng::seed_from_u64(seed);
let mut v: Vec<f64> = (0..n).map(|_| rng.random::<f64>() - 0.5).collect();
normalize_inplace(&mut v);
let mut lambda = 0.0f64;
for _ in 0..max_iter {
let mv = matvec(m, &v);
let new_lambda: f64 = mv.iter().zip(v.iter()).map(|(a, b)| a * b).sum();
let mut new_v = mv;
normalize_inplace(&mut new_v);
if (new_lambda - lambda).abs() < tol {
return (new_lambda, new_v);
}
lambda = new_lambda;
v = new_v;
}
(lambda, v)
}
fn inverse_iteration(
m: &Array2<f64>,
sigma: f64,
prev_evecs: &[Vec<f64>],
max_iter: usize,
tol: f64,
seed: u64,
) -> std::result::Result<(f64, Vec<f64>), String> {
let n = m.nrows();
if n == 0 {
return Err("empty matrix".into());
}
let mut a = m.to_owned();
for i in 0..n {
a[[i, i]] -= sigma;
}
let (lu, piv) = lu_decompose(&a)?;
let mut rng = StdRng::seed_from_u64(seed);
let mut v: Vec<f64> = (0..n).map(|_| rng.random::<f64>() - 0.5).collect();
deflate(&mut v, prev_evecs);
normalize_inplace(&mut v);
let mut eigenvalue = sigma;
for _ in 0..max_iter {
let mut w = lu_solve(&lu, &piv, &v)?;
deflate(&mut w, prev_evecs);
let norm = normalize_inplace(&mut w);
if norm < 1e-300 {
break;
}
let mv = matvec(m, &w);
let new_eigenvalue: f64 = mv.iter().zip(w.iter()).map(|(a, b)| a * b).sum();
let diff: f64 = w.iter().zip(v.iter()).map(|(a, b)| (a - b).powi(2)).sum::<f64>().sqrt();
v = w;
if (new_eigenvalue - eigenvalue).abs() < tol && diff < tol {
eigenvalue = new_eigenvalue;
break;
}
eigenvalue = new_eigenvalue;
}
Ok((eigenvalue, v))
}
fn lu_decompose(a: &Array2<f64>) -> std::result::Result<(Vec<Vec<f64>>, Vec<usize>), String> {
let n = a.nrows();
let mut lu: Vec<Vec<f64>> = (0..n).map(|i| a.row(i).to_vec()).collect();
let mut piv: Vec<usize> = (0..n).collect();
for k in 0..n {
let mut max_val = lu[k][k].abs();
let mut max_row = k;
for i in (k + 1)..n {
if lu[i][k].abs() > max_val {
max_val = lu[i][k].abs();
max_row = i;
}
}
if max_val < 1e-300 {
return Err(format!("singular matrix at column {k}"));
}
lu.swap(k, max_row);
piv.swap(k, max_row);
for i in (k + 1)..n {
lu[i][k] /= lu[k][k];
for j in (k + 1)..n {
let factor = lu[i][k] * lu[k][j];
lu[i][j] -= factor;
}
}
}
Ok((lu, piv))
}
fn lu_solve(
lu: &[Vec<f64>],
piv: &[usize],
b: &[f64],
) -> std::result::Result<Vec<f64>, String> {
let n = lu.len();
let mut x: Vec<f64> = vec![0.0; n];
for i in 0..n {
x[i] = b[piv[i]];
}
for i in 1..n {
for j in 0..i {
x[i] -= lu[i][j] * x[j];
}
}
for i in (0..n).rev() {
for j in (i + 1)..n {
x[i] -= lu[i][j] * x[j];
}
if lu[i][i].abs() < 1e-300 {
return Err(format!("zero pivot at {i}"));
}
x[i] /= lu[i][i];
}
Ok(x)
}
fn smallest_k_eigen(
m: &Array2<f64>,
k: usize,
seed: u64,
) -> std::result::Result<(Vec<f64>, Vec<Vec<f64>>), String> {
let n = m.nrows();
if k == 0 {
return Ok((vec![], vec![]));
}
if k > n {
return Err(format!("k={k} > n={n}"));
}
let (lambda_max, _) = power_iteration(m, 200, 1e-8, seed);
let mut eigenvalues: Vec<f64> = Vec::with_capacity(k);
let mut eigenvectors: Vec<Vec<f64>> = Vec::with_capacity(k);
let mut shift = -1e-4;
for idx in 0..k {
let trial_shift = if idx == 0 {
-1e-4 } else {
eigenvalues[idx - 1] - 1e-3
};
shift = trial_shift.min(shift);
let result = inverse_iteration(m, shift, &eigenvectors, 150, 1e-8, seed + idx as u64);
let (eval, evec) = match result {
Ok(r) => r,
Err(_) => {
inverse_iteration(
m,
shift + 1e-6,
&eigenvectors,
150,
1e-8,
seed + idx as u64 + 1000,
)?
}
};
let clamped = if eval.abs() < 1e-6 { 0.0 } else { eval };
eigenvalues.push(clamped);
eigenvectors.push(evec);
shift = clamped + (lambda_max - clamped) * 0.01;
}
Ok((eigenvalues, eigenvectors))
}
pub fn fiedler_vector(adj: &Array2<f64>) -> Result<(f64, Vec<f64>)> {
let n = adj.nrows();
if n < 2 {
return Err(GraphError::InvalidGraph(
"need at least 2 nodes for Fiedler vector".into(),
));
}
let lap = graph_laplacian(adj);
let (evals, evecs) = smallest_k_eigen(&lap, 2, 12345).map_err(|e| GraphError::LinAlgError {
operation: "fiedler_vector".into(),
details: e,
})?;
if evals.len() < 2 {
return Err(GraphError::LinAlgError {
operation: "fiedler_vector".into(),
details: "could not compute second eigenvalue".into(),
});
}
Ok((evals[1], evecs[1].clone()))
}
#[derive(Debug, Clone)]
pub struct SpectralClusterResult {
pub assignments: Vec<usize>,
pub eigenvalues: Vec<f64>,
pub embedding: Array2<f64>,
}
pub fn spectral_clustering(adj: &Array2<f64>, k: usize, seed: u64) -> Result<SpectralClusterResult> {
let n = adj.nrows();
if n == 0 {
return Err(GraphError::InvalidGraph("empty adjacency matrix".into()));
}
if k == 0 || k > n {
return Err(GraphError::InvalidParameter {
param: "k".into(),
value: k.to_string(),
expected: "1..=n".into(),
context: "spectral_clustering".into(),
});
}
let l_norm = normalized_laplacian(adj)?;
let (evals, evecs) =
smallest_k_eigen(&l_norm, k, seed).map_err(|e| GraphError::LinAlgError {
operation: "spectral_clustering".into(),
details: e,
})?;
let mut embedding = Array2::zeros((n, k));
for col in 0..k {
if col < evecs.len() {
for row in 0..n {
embedding[[row, col]] = evecs[col][row];
}
}
}
let assignments = kmeans_lloyd(&embedding, k, 300, seed)?;
Ok(SpectralClusterResult {
assignments,
eigenvalues: evals,
embedding,
})
}
fn kmeans_lloyd(data: &Array2<f64>, k: usize, max_iter: usize, seed: u64) -> Result<Vec<usize>> {
let n = data.nrows();
let d = data.ncols();
if n == 0 || k == 0 {
return Ok(vec![]);
}
let mut rng = StdRng::seed_from_u64(seed);
let mut centroids: Vec<Vec<f64>> = Vec::with_capacity(k);
let first_idx = rng.random_range(0..n);
centroids.push(data.row(first_idx).to_vec());
for _ in 1..k {
let dists: Vec<f64> = (0..n)
.map(|i| {
let row = data.row(i);
centroids
.iter()
.map(|c| {
row.iter()
.zip(c.iter())
.map(|(a, b)| (a - b).powi(2))
.sum::<f64>()
})
.fold(f64::INFINITY, f64::min)
})
.collect();
let total: f64 = dists.iter().sum();
let mut r = rng.random::<f64>() * total;
let mut chosen = n - 1;
for (i, &d_val) in dists.iter().enumerate() {
r -= d_val;
if r <= 0.0 {
chosen = i;
break;
}
}
centroids.push(data.row(chosen).to_vec());
}
let mut assignments = vec![0usize; n];
for _iter in 0..max_iter {
let mut changed = false;
for i in 0..n {
let row = data.row(i);
let mut best_dist = f64::INFINITY;
let mut best_c = 0;
for (ci, centroid) in centroids.iter().enumerate() {
let dist: f64 = row
.iter()
.zip(centroid.iter())
.map(|(a, b)| (a - b).powi(2))
.sum();
if dist < best_dist {
best_dist = dist;
best_c = ci;
}
}
if assignments[i] != best_c {
changed = true;
assignments[i] = best_c;
}
}
if !changed {
break;
}
let mut sums = vec![vec![0.0f64; d]; k];
let mut counts = vec![0usize; k];
for i in 0..n {
let c = assignments[i];
counts[c] += 1;
for j in 0..d {
sums[c][j] += data[[i, j]];
}
}
for c in 0..k {
if counts[c] > 0 {
for j in 0..d {
centroids[c][j] = sums[c][j] / counts[c] as f64;
}
}
}
}
Ok(assignments)
}
pub fn graph_fourier_transform(adj: &Array2<f64>, signal: &[f64]) -> Result<Vec<f64>> {
let n = adj.nrows();
if n == 0 {
return Err(GraphError::InvalidGraph("empty adjacency matrix".into()));
}
if signal.len() != n {
return Err(GraphError::InvalidParameter {
param: "signal".into(),
value: signal.len().to_string(),
expected: format!("{n}"),
context: "graph_fourier_transform".into(),
});
}
let lap = graph_laplacian(adj);
let (_, evecs) = smallest_k_eigen(&lap, n, 42).map_err(|e| GraphError::LinAlgError {
operation: "graph_fourier_transform".into(),
details: e,
})?;
let coeffs: Vec<f64> = evecs
.iter()
.map(|uk| uk.iter().zip(signal.iter()).map(|(a, b)| a * b).sum())
.collect();
Ok(coeffs)
}
fn pinv_symmetric(m: &Array2<f64>, tol: f64) -> std::result::Result<Array2<f64>, String> {
let n = m.nrows();
if n == 0 {
return Ok(Array2::zeros((0, 0)));
}
let (evals, evecs) = smallest_k_eigen(m, n, 999)?;
let mut pinv = Array2::zeros((n, n));
for (k, &lambda) in evals.iter().enumerate() {
if lambda.abs() <= tol {
continue;
}
let uk = &evecs[k];
for i in 0..n {
for j in 0..n {
pinv[[i, j]] += uk[i] * uk[j] / lambda;
}
}
}
Ok(pinv)
}
pub fn effective_resistance(adj: &Array2<f64>, i: usize, j: usize) -> Result<f64> {
let n = adj.nrows();
if n == 0 {
return Err(GraphError::InvalidGraph("empty adjacency matrix".into()));
}
if i >= n || j >= n {
return Err(GraphError::InvalidParameter {
param: "i or j".into(),
value: format!("({i},{j})"),
expected: format!("< {n}"),
context: "effective_resistance".into(),
});
}
if i == j {
return Ok(0.0);
}
let lap = graph_laplacian(adj);
let lp = pinv_symmetric(&lap, 1e-9).map_err(|e| GraphError::LinAlgError {
operation: "effective_resistance".into(),
details: e,
})?;
let r = lp[[i, i]] + lp[[j, j]] - 2.0 * lp[[i, j]];
Ok(r.max(0.0))
}
pub fn resistance_matrix(adj: &Array2<f64>) -> Result<Array2<f64>> {
let n = adj.nrows();
if n == 0 {
return Err(GraphError::InvalidGraph("empty adjacency matrix".into()));
}
if adj.ncols() != n {
return Err(GraphError::InvalidGraph("adjacency matrix must be square".into()));
}
let lap = graph_laplacian(adj);
let lp = pinv_symmetric(&lap, 1e-9).map_err(|e| GraphError::LinAlgError {
operation: "resistance_matrix".into(),
details: e,
})?;
let mut r_mat = Array2::zeros((n, n));
for i in 0..n {
for j in 0..n {
if i != j {
let r = lp[[i, i]] + lp[[j, j]] - 2.0 * lp[[i, j]];
r_mat[[i, j]] = r.max(0.0);
}
}
}
Ok(r_mat)
}
#[cfg(test)]
mod tests {
use super::*;
use scirs2_core::ndarray::Array2;
fn path_adj(n: usize) -> Array2<f64> {
let mut adj = Array2::zeros((n, n));
for i in 0..(n - 1) {
adj[[i, i + 1]] = 1.0;
adj[[i + 1, i]] = 1.0;
}
adj
}
fn complete_adj(n: usize) -> Array2<f64> {
let mut adj = Array2::ones((n, n));
for i in 0..n {
adj[[i, i]] = 0.0;
}
adj
}
fn complete_bipartite(m: usize, n: usize) -> Array2<f64> {
let total = m + n;
let mut adj = Array2::zeros((total, total));
for i in 0..m {
for j in m..total {
adj[[i, j]] = 1.0;
adj[[j, i]] = 1.0;
}
}
adj
}
#[test]
fn test_laplacian_path3() {
let adj = path_adj(3);
let l = graph_laplacian(&adj);
for i in 0..3 {
let row_sum: f64 = l.row(i).sum();
assert!(row_sum.abs() < 1e-12, "row {i} sum = {row_sum}");
}
assert!((l[[0, 0]] - 1.0).abs() < 1e-12);
assert!((l[[1, 1]] - 2.0).abs() < 1e-12);
assert!((l[[2, 2]] - 1.0).abs() < 1e-12);
}
#[test]
fn test_laplacian_complete4() {
let adj = complete_adj(4);
let l = graph_laplacian(&adj);
for i in 0..4 {
let row_sum: f64 = l.row(i).sum();
assert!(row_sum.abs() < 1e-12, "row {i} sum = {row_sum}");
assert!((l[[i, i]] - 3.0).abs() < 1e-12);
}
}
#[test]
fn test_normalized_laplacian_diagonal_ones() {
let adj = complete_adj(4); let l_norm = normalized_laplacian(&adj).expect("norm lap");
for i in 0..4 {
assert!(
(l_norm[[i, i]] - 1.0).abs() < 1e-10,
"diagonal[{i}] = {}",
l_norm[[i, i]]
);
}
}
#[test]
fn test_normalized_laplacian_symmetric() {
let adj = path_adj(5);
let l_norm = normalized_laplacian(&adj).expect("norm lap");
for i in 0..5 {
for j in 0..5 {
assert!(
(l_norm[[i, j]] - l_norm[[j, i]]).abs() < 1e-12,
"not symmetric at ({i},{j})"
);
}
}
}
#[test]
fn test_fiedler_value_positive_connected() {
let adj = path_adj(6);
let (lambda2, fv) = fiedler_vector(&adj).expect("fiedler");
assert!(lambda2 > 1e-6, "Fiedler value should be positive: {lambda2}");
assert_eq!(fv.len(), 6);
let norm: f64 = fv.iter().map(|x| x * x).sum::<f64>().sqrt();
assert!((norm - 1.0).abs() < 1e-6, "Fiedler vector should be unit norm: {norm}");
}
#[test]
fn test_fiedler_bisection() {
let adj = path_adj(4);
let (_lambda2, fv) = fiedler_vector(&adj).expect("fiedler");
let left_mean = (fv[0] + fv[1]) / 2.0;
let right_mean = (fv[2] + fv[3]) / 2.0;
assert!(
left_mean * right_mean < 0.0,
"Fiedler vector should bisect: left={left_mean}, right={right_mean}"
);
}
#[test]
fn test_fiedler_error_single_node() {
let adj = Array2::<f64>::zeros((1, 1));
assert!(fiedler_vector(&adj).is_err());
}
#[test]
fn test_spectral_clustering_shapes() {
let adj = path_adj(6);
let result = spectral_clustering(&adj, 2, 42).expect("spectral clustering");
assert_eq!(result.assignments.len(), 6);
assert_eq!(result.eigenvalues.len(), 2);
assert_eq!(result.embedding.nrows(), 6);
assert_eq!(result.embedding.ncols(), 2);
}
#[test]
fn test_spectral_clustering_two_components() {
let mut adj = Array2::zeros((6, 6));
adj[[0, 1]] = 1.0; adj[[1, 0]] = 1.0;
adj[[1, 2]] = 1.0; adj[[2, 1]] = 1.0;
adj[[3, 4]] = 1.0; adj[[4, 3]] = 1.0;
adj[[4, 5]] = 1.0; adj[[5, 4]] = 1.0;
let result = spectral_clustering(&adj, 2, 7).expect("spectral clustering");
let c0 = result.assignments[0];
let c3 = result.assignments[3];
assert_ne!(c0, c3, "disconnected components should be in different clusters");
assert_eq!(result.assignments[0], result.assignments[1]);
assert_eq!(result.assignments[1], result.assignments[2]);
assert_eq!(result.assignments[3], result.assignments[4]);
assert_eq!(result.assignments[4], result.assignments[5]);
}
#[test]
fn test_spectral_clustering_invalid_k() {
let adj = path_adj(4);
assert!(spectral_clustering(&adj, 0, 0).is_err());
assert!(spectral_clustering(&adj, 5, 0).is_err());
}
#[test]
fn test_gft_length() {
let adj = path_adj(4);
let signal = vec![1.0, 0.0, 0.0, 0.0];
let coeffs = graph_fourier_transform(&adj, &signal).expect("gft");
assert_eq!(coeffs.len(), 4);
}
#[test]
fn test_gft_constant_signal() {
let adj = complete_adj(4);
let signal = vec![1.0; 4];
let coeffs = graph_fourier_transform(&adj, &signal).expect("gft");
assert_eq!(coeffs.len(), 4);
let max_coeff = coeffs.iter().map(|x| x.abs()).fold(0.0f64, f64::max);
assert!(max_coeff > 0.5, "DC component should dominate: {coeffs:?}");
}
#[test]
fn test_effective_resistance_self() {
let adj = path_adj(4);
let r = effective_resistance(&adj, 1, 1).expect("eff res");
assert_eq!(r, 0.0);
}
#[test]
fn test_effective_resistance_path3() {
let adj = path_adj(3);
let r01 = effective_resistance(&adj, 0, 1).expect("r01");
let r12 = effective_resistance(&adj, 1, 2).expect("r12");
let r02 = effective_resistance(&adj, 0, 2).expect("r02");
assert!((r01 - 1.0).abs() < 1e-4, "R(0,1) = {r01}");
assert!((r12 - 1.0).abs() < 1e-4, "R(1,2) = {r12}");
assert!((r02 - 2.0).abs() < 1e-4, "R(0,2) = {r02}");
}
#[test]
fn test_effective_resistance_complete_graph() {
let n = 4;
let adj = complete_adj(n);
let r = effective_resistance(&adj, 0, 1).expect("eff res");
let expected = 2.0 / n as f64;
assert!(
(r - expected).abs() < 1e-4,
"K{n} effective resistance = {r}, expected {expected}"
);
}
#[test]
fn test_resistance_matrix_shape() {
let adj = path_adj(4);
let r_mat = resistance_matrix(&adj).expect("res mat");
assert_eq!(r_mat.nrows(), 4);
assert_eq!(r_mat.ncols(), 4);
}
#[test]
fn test_resistance_matrix_symmetric() {
let adj = path_adj(5);
let r_mat = resistance_matrix(&adj).expect("res mat");
for i in 0..5 {
for j in 0..5 {
assert!(
(r_mat[[i, j]] - r_mat[[j, i]]).abs() < 1e-8,
"not symmetric at ({i},{j})"
);
}
}
}
#[test]
fn test_resistance_matrix_zero_diagonal() {
let adj = complete_adj(4);
let r_mat = resistance_matrix(&adj).expect("res mat");
for i in 0..4 {
assert!(r_mat[[i, i]].abs() < 1e-10, "diagonal[{i}] = {}", r_mat[[i, i]]);
}
}
#[test]
fn test_resistance_matrix_path3_values() {
let adj = path_adj(3);
let r_mat = resistance_matrix(&adj).expect("res mat");
assert!((r_mat[[0, 1]] - 1.0).abs() < 1e-4, "R(0,1) = {}", r_mat[[0, 1]]);
assert!((r_mat[[1, 2]] - 1.0).abs() < 1e-4, "R(1,2) = {}", r_mat[[1, 2]]);
assert!((r_mat[[0, 2]] - 2.0).abs() < 1e-4, "R(0,2) = {}", r_mat[[0, 2]]);
}
#[test]
fn test_bipartite_laplacian_max_eigenvalue() {
let adj = complete_bipartite(2, 2);
let l_norm = normalized_laplacian(&adj).expect("norm lap");
let (lambda_max, _) = power_iteration(&l_norm, 500, 1e-10, 77);
assert!(
(lambda_max - 2.0).abs() < 0.05,
"K_{{2,2}} max eigenvalue of L_norm = {lambda_max}, expected 2.0"
);
}
}