use crate::array::Array;
use crate::error::{NumRs2Error, Result};
use num_traits::Float;
use std::fmt::Debug;
#[derive(Debug, Clone)]
pub struct TuckerResult<T: Clone> {
pub core: Array<T>,
pub factors: Vec<Array<T>>,
pub original_shape: Vec<usize>,
pub ranks: Vec<usize>,
}
#[derive(Debug, Clone)]
pub struct CpResult<T: Clone> {
pub factors: Vec<Array<T>>,
pub weights: Vec<T>,
pub rank: usize,
pub error: T,
pub iterations: usize,
}
#[derive(Debug, Clone)]
pub struct DecompConfig {
pub max_iterations: usize,
pub tolerance: f64,
pub normalize: bool,
}
impl Default for DecompConfig {
fn default() -> Self {
DecompConfig {
max_iterations: 100,
tolerance: 1e-6,
normalize: true,
}
}
}
pub fn tucker_decomposition<T>(tensor: &Array<T>, ranks: &[usize]) -> Result<TuckerResult<T>>
where
T: Float + Clone + Debug + Default + Send + Sync + std::iter::Sum,
{
let shape = tensor.shape();
let ndim = shape.len();
if ranks.len() != ndim {
return Err(NumRs2Error::DimensionMismatch(format!(
"Number of ranks ({}) must match tensor dimensions ({})",
ranks.len(),
ndim
)));
}
for (i, (&rank, &dim)) in ranks.iter().zip(shape.iter()).enumerate() {
if rank > dim {
return Err(NumRs2Error::ValueError(format!(
"Rank {} for mode {} exceeds dimension {}",
rank, i, dim
)));
}
if rank == 0 {
return Err(NumRs2Error::ValueError(format!(
"Rank for mode {} cannot be zero",
i
)));
}
}
let mut factors = Vec::with_capacity(ndim);
for mode in 0..ndim {
let unfolded = mode_unfold(tensor, mode)?;
let (u, _, _) = truncated_svd(&unfolded, ranks[mode])?;
factors.push(u);
}
let mut core = tensor.clone();
for (mode, factor) in factors.iter().enumerate() {
let factor_t = transpose_2d(factor)?;
core = n_mode_product(&core, &factor_t, mode)?;
}
Ok(TuckerResult {
core,
factors,
original_shape: shape.clone(),
ranks: ranks.to_vec(),
})
}
pub fn tucker_reconstruct<T>(result: &TuckerResult<T>) -> Result<Array<T>>
where
T: Float + Clone + Debug + Default + Send + Sync + std::iter::Sum,
{
let mut reconstructed = result.core.clone();
for (mode, factor) in result.factors.iter().enumerate() {
reconstructed = n_mode_product(&reconstructed, factor, mode)?;
}
Ok(reconstructed)
}
pub fn cp_als<T>(
tensor: &Array<T>,
rank: usize,
max_iter: usize,
tolerance: T,
) -> Result<CpResult<T>>
where
T: Float + Clone + Debug + Default + Send + Sync + std::iter::Sum,
{
let shape = tensor.shape();
let ndim = shape.len();
if rank == 0 {
return Err(NumRs2Error::ValueError("Rank must be positive".to_string()));
}
let mut factors: Vec<Array<T>> = shape
.iter()
.map(|&dim| initialize_factor(dim, rank))
.collect();
let mut prev_error = T::infinity();
for iteration in 0..max_iter {
for mode in 0..ndim {
let new_factor = update_factor_als(tensor, &factors, mode)?;
factors[mode] = new_factor;
}
let reconstructed = cp_reconstruct_from_factors(&factors)?;
let error = frobenius_error(tensor, &reconstructed)?;
let error_change = (prev_error - error).abs();
if error_change < tolerance {
let weights = extract_weights(&factors);
return Ok(CpResult {
factors,
weights,
rank,
error,
iterations: iteration + 1,
});
}
prev_error = error;
}
let weights = extract_weights(&factors);
let error = prev_error;
Ok(CpResult {
factors,
weights,
rank,
error,
iterations: max_iter,
})
}
pub fn cp_als_with_config<T>(
tensor: &Array<T>,
rank: usize,
config: &DecompConfig,
) -> Result<CpResult<T>>
where
T: Float + Clone + Debug + Default + Send + Sync + std::iter::Sum,
{
cp_als(
tensor,
rank,
config.max_iterations,
T::from(config.tolerance).unwrap_or(T::epsilon()),
)
}
pub fn cp_reconstruct<T>(result: &CpResult<T>) -> Result<Array<T>>
where
T: Float + Clone + Debug + Default + Send + Sync + std::iter::Sum,
{
cp_reconstruct_from_factors(&result.factors)
}
pub fn nonnegative_cp_als<T>(
tensor: &Array<T>,
rank: usize,
max_iter: usize,
tolerance: T,
) -> Result<CpResult<T>>
where
T: Float + Clone + Debug + Default + Send + Sync + std::iter::Sum,
{
let shape = tensor.shape();
let ndim = shape.len();
for val in tensor.to_vec() {
if val < T::zero() {
return Err(NumRs2Error::ValueError(
"Tensor must be non-negative for non-negative CP decomposition".to_string(),
));
}
}
let mut factors: Vec<Array<T>> = shape
.iter()
.map(|&dim| initialize_nonnegative_factor(dim, rank))
.collect();
let mut prev_error = T::infinity();
for iteration in 0..max_iter {
for mode in 0..ndim {
let new_factor = update_factor_als_nonnegative(tensor, &factors, mode)?;
factors[mode] = new_factor;
}
let reconstructed = cp_reconstruct_from_factors(&factors)?;
let error = frobenius_error(tensor, &reconstructed)?;
let error_change = (prev_error - error).abs();
if error_change < tolerance {
let weights = extract_weights(&factors);
return Ok(CpResult {
factors,
weights,
rank,
error,
iterations: iteration + 1,
});
}
prev_error = error;
}
let weights = extract_weights(&factors);
Ok(CpResult {
factors,
weights,
rank,
error: prev_error,
iterations: max_iter,
})
}
fn mode_unfold<T>(tensor: &Array<T>, mode: usize) -> Result<Array<T>>
where
T: Float + Clone + Debug + Default + Send + Sync,
{
use scirs2_core::parallel_ops::*;
const PARALLEL_THRESHOLD: usize = 1000;
let shape = tensor.shape();
let ndim = shape.len();
if mode >= ndim {
return Err(NumRs2Error::ValueError(format!(
"Mode {} exceeds tensor dimensions {}",
mode, ndim
)));
}
let mode_size = shape[mode];
let other_size: usize = shape
.iter()
.enumerate()
.filter(|(i, _)| *i != mode)
.map(|(_, &s)| s)
.product();
let data = tensor.to_vec();
let mut strides = vec![1usize; ndim];
for i in (0..ndim - 1).rev() {
strides[i] = strides[i + 1] * shape[i + 1];
}
let mut col_strides = Vec::with_capacity(ndim);
let mut cumulative = 1usize;
for d in (0..ndim).rev() {
if d != mode {
col_strides.push((d, cumulative));
cumulative *= shape[d];
}
}
col_strides.reverse();
if is_parallel_enabled() && data.len() >= PARALLEL_THRESHOLD {
let results: Vec<(usize, T)> = (0..data.len())
.into_par_iter()
.map(|idx| {
let row = (idx / strides[mode]) % shape[mode];
let mut col = 0usize;
for &(d, col_stride) in &col_strides {
let dim_idx = (idx / strides[d]) % shape[d];
col += dim_idx * col_stride;
}
let result_idx = row * other_size + col;
(result_idx, data[idx])
})
.collect();
let mut result = vec![T::zero(); mode_size * other_size];
for (idx, val) in results {
result[idx] = val;
}
Ok(Array::from_vec(result).reshape(&[mode_size, other_size]))
} else {
let mut result = vec![T::zero(); mode_size * other_size];
for idx in 0..data.len() {
let row = (idx / strides[mode]) % shape[mode];
let mut col = 0usize;
for &(d, col_stride) in &col_strides {
let dim_idx = (idx / strides[d]) % shape[d];
col += dim_idx * col_stride;
}
result[row * other_size + col] = data[idx];
}
Ok(Array::from_vec(result).reshape(&[mode_size, other_size]))
}
}
fn mode_fold<T>(matrix: &Array<T>, mode: usize, shape: &[usize]) -> Result<Array<T>>
where
T: Float + Clone + Debug + Default + Send + Sync,
{
use scirs2_core::parallel_ops::*;
const PARALLEL_THRESHOLD: usize = 1000;
let ndim = shape.len();
let total_size: usize = shape.iter().product();
let mat_shape = matrix.shape();
if mat_shape.len() != 2 {
return Err(NumRs2Error::DimensionMismatch(
"Input must be a 2D matrix".to_string(),
));
}
let mode_size = mat_shape[0];
let other_size = mat_shape[1];
if mode_size != shape[mode] {
return Err(NumRs2Error::DimensionMismatch(format!(
"Matrix row count {} doesn't match mode dimension {}",
mode_size, shape[mode]
)));
}
let mat_data = matrix.to_vec();
let mut strides = vec![1usize; ndim];
for i in (0..ndim - 1).rev() {
strides[i] = strides[i + 1] * shape[i + 1];
}
let mut col_divisors_and_strides = Vec::with_capacity(ndim);
let mut col_stride = other_size;
for d in (0..ndim).rev() {
if d != mode {
col_stride /= shape[d];
col_divisors_and_strides.push((col_stride, shape[d], strides[d]));
}
}
col_divisors_and_strides.reverse();
let mode_stride = strides[mode];
let total_elements = mode_size * other_size;
if is_parallel_enabled() && total_elements >= PARALLEL_THRESHOLD {
let results: Vec<(usize, T)> = (0..total_elements)
.into_par_iter()
.map(|idx| {
let row = idx / other_size;
let col = idx % other_size;
let mut flat_idx = row * mode_stride;
let mut remaining = col;
for &(divisor, dim_size, stride) in &col_divisors_and_strides {
let dim_idx = (remaining / divisor) % dim_size;
flat_idx += dim_idx * stride;
remaining %= divisor;
}
(flat_idx, mat_data[row * other_size + col])
})
.collect();
let mut result = vec![T::zero(); total_size];
for (idx, val) in results {
result[idx] = val;
}
Ok(Array::from_vec(result).reshape(shape))
} else {
let mut result = vec![T::zero(); total_size];
for row in 0..mode_size {
let row_contrib = row * mode_stride;
for col in 0..other_size {
let mut flat_idx = row_contrib;
let mut remaining = col;
for &(divisor, dim_size, stride) in &col_divisors_and_strides {
let dim_idx = (remaining / divisor) % dim_size;
flat_idx += dim_idx * stride;
remaining %= divisor;
}
result[flat_idx] = mat_data[row * other_size + col];
}
}
Ok(Array::from_vec(result).reshape(shape))
}
}
fn n_mode_product<T>(tensor: &Array<T>, matrix: &Array<T>, mode: usize) -> Result<Array<T>>
where
T: Float + Clone + Debug + Default + Send + Sync + std::iter::Sum,
{
let t_shape = tensor.shape();
let m_shape = matrix.shape();
if m_shape.len() != 2 {
return Err(NumRs2Error::DimensionMismatch(
"Matrix must be 2D".to_string(),
));
}
let unfolded = mode_unfold(tensor, mode)?;
let product = matrix_multiply(matrix, &unfolded)?;
let mut new_shape = t_shape.clone();
new_shape[mode] = m_shape[0];
mode_fold(&product, mode, &new_shape)
}
fn truncated_svd<T>(matrix: &Array<T>, k: usize) -> Result<(Array<T>, Array<T>, Array<T>)>
where
T: Float + Clone + Debug + Default + Send + Sync + std::iter::Sum,
{
let shape = matrix.shape();
if shape.len() != 2 {
return Err(NumRs2Error::DimensionMismatch(
"SVD requires a 2D matrix".to_string(),
));
}
let (m, n) = (shape[0], shape[1]);
let rank = k.min(m).min(n);
let mut v = gram_schmidt(&random_matrix(n, rank))?;
let max_iter = 50;
for _ in 0..max_iter {
let av = matrix_multiply(matrix, &v)?;
let (u, _) = qr_decomposition(&av)?;
let at = transpose_2d(matrix)?;
let atu = matrix_multiply(&at, &u)?;
let (v_new, _) = qr_decomposition(&atu)?;
v = v_new;
}
let av = matrix_multiply(matrix, &v)?;
let (u, r) = qr_decomposition(&av)?;
let mut s_data = vec![T::zero(); rank];
for i in 0..rank {
s_data[i] = r.get(&[i, i]).unwrap_or(T::zero()).abs();
}
let vt = transpose_2d(&v)?;
Ok((u, Array::from_vec(s_data), vt))
}
fn transpose_2d<T>(matrix: &Array<T>) -> Result<Array<T>>
where
T: Float + Clone + Debug + Default,
{
let shape = matrix.shape();
if shape.len() != 2 {
return Err(NumRs2Error::DimensionMismatch(
"Transpose requires a 2D matrix".to_string(),
));
}
let (m, n) = (shape[0], shape[1]);
let data = matrix.to_vec();
let mut result = vec![T::zero(); m * n];
for i in 0..m {
for j in 0..n {
result[j * m + i] = data[i * n + j];
}
}
Ok(Array::from_vec(result).reshape(&[n, m]))
}
fn matrix_multiply<T>(a: &Array<T>, b: &Array<T>) -> Result<Array<T>>
where
T: Float + Clone + Debug + Default + std::iter::Sum,
{
let a_shape = a.shape();
let b_shape = b.shape();
if a_shape.len() != 2 || b_shape.len() != 2 {
return Err(NumRs2Error::DimensionMismatch(
"Matrix multiplication requires 2D matrices".to_string(),
));
}
if a_shape[1] != b_shape[0] {
return Err(NumRs2Error::ShapeMismatch {
expected: vec![a_shape[0], b_shape[0]],
actual: vec![a_shape[1]],
});
}
let (m, k) = (a_shape[0], a_shape[1]);
let n = b_shape[1];
let a_data = a.to_vec();
let b_data = b.to_vec();
let mut result = vec![T::zero(); m * n];
for i in 0..m {
for j in 0..n {
let mut sum = T::zero();
for l in 0..k {
sum = sum + a_data[i * k + l] * b_data[l * n + j];
}
result[i * n + j] = sum;
}
}
Ok(Array::from_vec(result).reshape(&[m, n]))
}
fn initialize_factor<T>(rows: usize, cols: usize) -> Array<T>
where
T: Float + Default,
{
let mut data = vec![T::zero(); rows * cols];
let mut seed = 42u64;
for val in data.iter_mut() {
seed = seed
.wrapping_mul(6364136223846793005)
.wrapping_add(1442695040888963407);
let rand = ((seed >> 33) as f64) / (u32::MAX as f64);
*val =
T::from(rand).unwrap_or(T::one() / T::from(2.0).expect("2.0 is a valid f64 constant"));
}
Array::from_vec(data).reshape(&[rows, cols])
}
fn initialize_nonnegative_factor<T>(rows: usize, cols: usize) -> Array<T>
where
T: Float + Default,
{
let mut data = vec![T::zero(); rows * cols];
let mut seed = 42u64;
for val in data.iter_mut() {
seed = seed
.wrapping_mul(6364136223846793005)
.wrapping_add(1442695040888963407);
let rand = ((seed >> 33) as f64) / (u32::MAX as f64);
*val = T::from(rand.abs())
.unwrap_or(T::one() / T::from(2.0).expect("2.0 is a valid f64 constant"));
}
Array::from_vec(data).reshape(&[rows, cols])
}
fn random_matrix<T>(rows: usize, cols: usize) -> Array<T>
where
T: Float + Default,
{
initialize_factor(rows, cols)
}
fn gram_schmidt<T>(matrix: &Array<T>) -> Result<Array<T>>
where
T: Float + Clone + Debug + Default + std::iter::Sum,
{
let shape = matrix.shape();
if shape.len() != 2 {
return Err(NumRs2Error::DimensionMismatch(
"Gram-Schmidt requires 2D matrix".to_string(),
));
}
let (m, n) = (shape[0], shape[1]);
let data = matrix.to_vec();
let mut result = vec![T::zero(); m * n];
for j in 0..n {
for i in 0..m {
result[i * n + j] = data[i * n + j];
}
for k in 0..j {
let mut dot = T::zero();
let mut norm_sq = T::zero();
for i in 0..m {
dot = dot + result[i * n + k] * data[i * n + j];
norm_sq = norm_sq + result[i * n + k] * result[i * n + k];
}
if norm_sq > T::epsilon() {
let coef = dot / norm_sq;
for i in 0..m {
result[i * n + j] = result[i * n + j] - coef * result[i * n + k];
}
}
}
let mut norm_sq = T::zero();
for i in 0..m {
norm_sq = norm_sq + result[i * n + j] * result[i * n + j];
}
let norm = norm_sq.sqrt();
if norm > T::epsilon() {
for i in 0..m {
result[i * n + j] = result[i * n + j] / norm;
}
}
}
Ok(Array::from_vec(result).reshape(&[m, n]))
}
fn qr_decomposition<T>(matrix: &Array<T>) -> Result<(Array<T>, Array<T>)>
where
T: Float + Clone + Debug + Default + std::iter::Sum,
{
let shape = matrix.shape();
if shape.len() != 2 {
return Err(NumRs2Error::DimensionMismatch(
"QR requires 2D matrix".to_string(),
));
}
let (m, n) = (shape[0], shape[1]);
let q = gram_schmidt(matrix)?;
let q_t = transpose_2d(&q)?;
let r = matrix_multiply(&q_t, matrix)?;
Ok((q, r))
}
fn update_factor_als<T>(tensor: &Array<T>, factors: &[Array<T>], mode: usize) -> Result<Array<T>>
where
T: Float + Clone + Debug + Default + Send + Sync + std::iter::Sum,
{
let shape = tensor.shape();
let ndim = shape.len();
let rank = factors[0].shape()[1];
let unfolded = mode_unfold(tensor, mode)?;
let kr = khatri_rao_except(factors, mode)?;
let kr_t = transpose_2d(&kr)?;
let gram = matrix_multiply(&kr_t, &kr)?;
let gram_reg = add_regularization(
&gram,
T::from(1e-10).expect("1e-10 is a valid f64 constant"),
)?;
let gram_inv = inverse_2x2_or_general(&gram_reg)?;
let temp = matrix_multiply(&unfolded, &kr)?;
let new_factor = matrix_multiply(&temp, &gram_inv)?;
normalize_columns(&new_factor)
}
fn update_factor_als_nonnegative<T>(
tensor: &Array<T>,
factors: &[Array<T>],
mode: usize,
) -> Result<Array<T>>
where
T: Float + Clone + Debug + Default + Send + Sync + std::iter::Sum,
{
let factor = update_factor_als(tensor, factors, mode)?;
let data = factor.to_vec();
let projected: Vec<T> = data.iter().map(|&x| x.max(T::zero())).collect();
Ok(Array::from_vec(projected).reshape(&factor.shape()))
}
fn khatri_rao<T>(a: &Array<T>, b: &Array<T>) -> Result<Array<T>>
where
T: Float + Clone + Debug + Default,
{
let a_shape = a.shape();
let b_shape = b.shape();
if a_shape.len() != 2 || b_shape.len() != 2 {
return Err(NumRs2Error::DimensionMismatch(
"Khatri-Rao requires 2D matrices".to_string(),
));
}
if a_shape[1] != b_shape[1] {
return Err(NumRs2Error::DimensionMismatch(
"Khatri-Rao requires same number of columns".to_string(),
));
}
let (m, r) = (a_shape[0], a_shape[1]);
let n = b_shape[0];
let a_data = a.to_vec();
let b_data = b.to_vec();
let mut result = vec![T::zero(); m * n * r];
for k in 0..r {
for i in 0..m {
for j in 0..n {
result[(i * n + j) * r + k] = a_data[i * r + k] * b_data[j * r + k];
}
}
}
Ok(Array::from_vec(result).reshape(&[m * n, r]))
}
fn khatri_rao_except<T>(factors: &[Array<T>], exclude: usize) -> Result<Array<T>>
where
T: Float + Clone + Debug + Default,
{
let ndim = factors.len();
let mut result: Option<Array<T>> = None;
for i in (0..ndim).rev() {
if i == exclude {
continue;
}
match result {
None => result = Some(factors[i].clone()),
Some(ref prev) => {
result = Some(khatri_rao(&factors[i], prev)?);
}
}
}
result.ok_or_else(|| {
NumRs2Error::ValueError("Need at least 2 factors for Khatri-Rao".to_string())
})
}
fn add_regularization<T>(matrix: &Array<T>, lambda: T) -> Result<Array<T>>
where
T: Float + Clone + Debug + Default,
{
let shape = matrix.shape();
if shape.len() != 2 || shape[0] != shape[1] {
return Err(NumRs2Error::DimensionMismatch(
"Regularization requires square matrix".to_string(),
));
}
let n = shape[0];
let mut data = matrix.to_vec();
for i in 0..n {
data[i * n + i] = data[i * n + i] + lambda;
}
Ok(Array::from_vec(data).reshape(&[n, n]))
}
fn inverse_2x2_or_general<T>(matrix: &Array<T>) -> Result<Array<T>>
where
T: Float + Clone + Debug + Default + std::iter::Sum,
{
let shape = matrix.shape();
if shape.len() != 2 || shape[0] != shape[1] {
return Err(NumRs2Error::DimensionMismatch(
"Inverse requires square matrix".to_string(),
));
}
let n = shape[0];
let data = matrix.to_vec();
if n == 1 {
if data[0].abs() < T::epsilon() {
return Err(NumRs2Error::InvalidOperation("Singular matrix".to_string()));
}
return Ok(Array::from_vec(vec![T::one() / data[0]]).reshape(&[1, 1]));
}
if n == 2 {
let det = data[0] * data[3] - data[1] * data[2];
if det.abs() < T::epsilon() {
return Err(NumRs2Error::InvalidOperation("Singular matrix".to_string()));
}
let inv_det = T::one() / det;
return Ok(Array::from_vec(vec![
data[3] * inv_det,
-data[1] * inv_det,
-data[2] * inv_det,
data[0] * inv_det,
])
.reshape(&[2, 2]));
}
gauss_jordan_inverse(matrix)
}
fn gauss_jordan_inverse<T>(matrix: &Array<T>) -> Result<Array<T>>
where
T: Float + Clone + Debug + Default,
{
let shape = matrix.shape();
let n = shape[0];
let data = matrix.to_vec();
let mut aug = vec![T::zero(); n * 2 * n];
for i in 0..n {
for j in 0..n {
aug[i * 2 * n + j] = data[i * n + j];
}
aug[i * 2 * n + n + i] = T::one();
}
for i in 0..n {
let mut max_row = i;
let mut max_val = aug[i * 2 * n + i].abs();
for k in i + 1..n {
let val = aug[k * 2 * n + i].abs();
if val > max_val {
max_val = val;
max_row = k;
}
}
if max_val < T::epsilon() {
return Err(NumRs2Error::InvalidOperation("Singular matrix".to_string()));
}
if max_row != i {
for j in 0..2 * n {
aug.swap(i * 2 * n + j, max_row * 2 * n + j);
}
}
let pivot = aug[i * 2 * n + i];
for j in 0..2 * n {
aug[i * 2 * n + j] = aug[i * 2 * n + j] / pivot;
}
for k in 0..n {
if k != i {
let factor = aug[k * 2 * n + i];
for j in 0..2 * n {
aug[k * 2 * n + j] = aug[k * 2 * n + j] - factor * aug[i * 2 * n + j];
}
}
}
}
let mut result = vec![T::zero(); n * n];
for i in 0..n {
for j in 0..n {
result[i * n + j] = aug[i * 2 * n + n + j];
}
}
Ok(Array::from_vec(result).reshape(&[n, n]))
}
fn normalize_columns<T>(matrix: &Array<T>) -> Result<Array<T>>
where
T: Float + Clone + Debug + Default + std::iter::Sum,
{
let shape = matrix.shape();
if shape.len() != 2 {
return Err(NumRs2Error::DimensionMismatch(
"Normalize requires 2D matrix".to_string(),
));
}
let (m, n) = (shape[0], shape[1]);
let mut data = matrix.to_vec();
for j in 0..n {
let mut norm_sq = T::zero();
for i in 0..m {
norm_sq = norm_sq + data[i * n + j] * data[i * n + j];
}
let norm = norm_sq.sqrt();
if norm > T::epsilon() {
for i in 0..m {
data[i * n + j] = data[i * n + j] / norm;
}
}
}
Ok(Array::from_vec(data).reshape(&[m, n]))
}
fn extract_weights<T>(factors: &[Array<T>]) -> Vec<T>
where
T: Float + Clone + Debug + Default + std::iter::Sum,
{
if factors.is_empty() {
return Vec::new();
}
let shape = factors[0].shape();
let rank = shape[1];
let mut weights = vec![T::one(); rank];
for factor in factors {
let data = factor.to_vec();
let f_shape = factor.shape();
let (m, n) = (f_shape[0], f_shape[1]);
for j in 0..n.min(rank) {
let mut norm_sq = T::zero();
for i in 0..m {
norm_sq = norm_sq + data[i * n + j] * data[i * n + j];
}
weights[j] = weights[j] * norm_sq.sqrt();
}
}
weights
}
fn cp_reconstruct_from_factors<T>(factors: &[Array<T>]) -> Result<Array<T>>
where
T: Float + Clone + Debug + Default + std::iter::Sum,
{
if factors.is_empty() {
return Err(NumRs2Error::ValueError(
"Need at least one factor".to_string(),
));
}
let ndim = factors.len();
let rank = factors[0].shape()[1];
let shape: Vec<usize> = factors.iter().map(|f| f.shape()[0]).collect();
let total_size: usize = shape.iter().product();
let mut strides = vec![1usize; ndim];
for i in (0..ndim - 1).rev() {
strides[i] = strides[i + 1] * shape[i + 1];
}
let mut result = vec![T::zero(); total_size];
let factor_data: Vec<Vec<T>> = factors.iter().map(|f| f.to_vec()).collect();
for r in 0..rank {
for idx in 0..total_size {
let mut multi_idx = vec![0usize; ndim];
let mut remaining = idx;
for d in 0..ndim {
multi_idx[d] = remaining / strides[d];
remaining %= strides[d];
}
let mut prod = T::one();
for (d, &m_idx) in multi_idx.iter().enumerate() {
let f_cols = factors[d].shape()[1];
prod = prod * factor_data[d][m_idx * f_cols + r];
}
result[idx] = result[idx] + prod;
}
}
Ok(Array::from_vec(result).reshape(&shape))
}
fn frobenius_error<T>(a: &Array<T>, b: &Array<T>) -> Result<T>
where
T: Float + Clone + Debug + Default + std::iter::Sum,
{
let a_data = a.to_vec();
let b_data = b.to_vec();
if a_data.len() != b_data.len() {
return Err(NumRs2Error::ShapeMismatch {
expected: a.shape(),
actual: b.shape(),
});
}
let sum: T = a_data
.iter()
.zip(b_data.iter())
.map(|(&x, &y)| (x - y) * (x - y))
.sum();
Ok(sum.sqrt())
}
#[cfg(test)]
mod tests {
use super::*;
#[test]
fn test_mode_unfold_2d() {
let tensor = Array::from_vec(vec![1.0, 2.0, 3.0, 4.0, 5.0, 6.0]).reshape(&[2, 3]);
let unfolded_0 = mode_unfold(&tensor, 0).expect("mode unfold should succeed");
assert_eq!(unfolded_0.shape(), vec![2, 3]);
let unfolded_1 = mode_unfold(&tensor, 1).expect("mode unfold should succeed");
assert_eq!(unfolded_1.shape(), vec![3, 2]);
}
#[test]
fn test_mode_unfold_3d() {
let tensor =
Array::from_vec(vec![1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0]).reshape(&[2, 2, 2]);
let unfolded_0 = mode_unfold(&tensor, 0).expect("mode unfold should succeed");
assert_eq!(unfolded_0.shape(), vec![2, 4]);
let unfolded_1 = mode_unfold(&tensor, 1).expect("mode unfold should succeed");
assert_eq!(unfolded_1.shape(), vec![2, 4]);
let unfolded_2 = mode_unfold(&tensor, 2).expect("mode unfold should succeed");
assert_eq!(unfolded_2.shape(), vec![2, 4]);
}
#[test]
fn test_transpose_2d() {
let matrix = Array::from_vec(vec![1.0, 2.0, 3.0, 4.0, 5.0, 6.0]).reshape(&[2, 3]);
let transposed = transpose_2d(&matrix).expect("transpose should succeed");
assert_eq!(transposed.shape(), vec![3, 2]);
assert_eq!(transposed.get(&[0, 0]).expect("valid index"), 1.0);
assert_eq!(transposed.get(&[0, 1]).expect("valid index"), 4.0);
assert_eq!(transposed.get(&[1, 0]).expect("valid index"), 2.0);
}
#[test]
fn test_matrix_multiply() {
let a = Array::from_vec(vec![1.0, 2.0, 3.0, 4.0]).reshape(&[2, 2]);
let b = Array::from_vec(vec![5.0, 6.0, 7.0, 8.0]).reshape(&[2, 2]);
let c = matrix_multiply(&a, &b).expect("matrix multiply should succeed");
assert_eq!(c.shape(), vec![2, 2]);
assert_eq!(c.get(&[0, 0]).expect("valid index"), 19.0); assert_eq!(c.get(&[0, 1]).expect("valid index"), 22.0); assert_eq!(c.get(&[1, 0]).expect("valid index"), 43.0); assert_eq!(c.get(&[1, 1]).expect("valid index"), 50.0); }
#[test]
fn test_tucker_decomposition_simple() {
let tensor =
Array::from_vec(vec![1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0]).reshape(&[2, 2, 2]);
let result =
tucker_decomposition(&tensor, &[2, 2, 2]).expect("tucker decomposition should succeed");
assert_eq!(result.core.shape(), vec![2, 2, 2]);
assert_eq!(result.factors.len(), 3);
assert_eq!(result.factors[0].shape(), vec![2, 2]);
assert_eq!(result.factors[1].shape(), vec![2, 2]);
assert_eq!(result.factors[2].shape(), vec![2, 2]);
}
#[test]
fn test_tucker_decomposition_reduced_rank() {
let tensor =
Array::from_vec(vec![1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0]).reshape(&[2, 2, 2]);
let result =
tucker_decomposition(&tensor, &[1, 1, 1]).expect("tucker decomposition should succeed");
assert_eq!(result.core.shape(), vec![1, 1, 1]);
assert_eq!(result.factors[0].shape(), vec![2, 1]);
assert_eq!(result.factors[1].shape(), vec![2, 1]);
assert_eq!(result.factors[2].shape(), vec![2, 1]);
}
#[test]
fn test_cp_als_simple() {
let tensor =
Array::from_vec(vec![1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0]).reshape(&[2, 2, 2]);
let result = cp_als(&tensor, 2, 50, 1e-4).expect("CP-ALS should succeed");
assert_eq!(result.rank, 2);
assert_eq!(result.factors.len(), 3);
assert_eq!(result.factors[0].shape(), vec![2, 2]);
assert_eq!(result.factors[1].shape(), vec![2, 2]);
assert_eq!(result.factors[2].shape(), vec![2, 2]);
}
#[test]
fn test_cp_reconstruct() {
let tensor =
Array::from_vec(vec![1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0]).reshape(&[2, 2, 2]);
let result = cp_als(&tensor, 3, 100, 1e-6).expect("CP-ALS should succeed");
let reconstructed = cp_reconstruct(&result).expect("CP reconstruct should succeed");
assert_eq!(reconstructed.shape(), tensor.shape());
let frobenius_norm: f64 = tensor.to_vec().iter().map(|x| x * x).sum::<f64>().sqrt();
let relative_error = result.error / frobenius_norm;
assert!(
relative_error < 2.0, "Relative error {} indicates algorithm failure",
relative_error
);
}
#[test]
fn test_khatri_rao() {
let a = Array::from_vec(vec![1.0, 2.0, 3.0, 4.0]).reshape(&[2, 2]);
let b = Array::from_vec(vec![5.0, 6.0, 7.0, 8.0, 9.0, 10.0]).reshape(&[3, 2]);
let kr = khatri_rao(&a, &b).expect("Khatri-Rao product should succeed");
assert_eq!(kr.shape(), vec![6, 2]);
}
#[test]
fn test_gram_schmidt() {
let matrix =
Array::from_vec(vec![1.0, 1.0, 1.0, 0.0, 1.0, 1.0, 1.0, 1.0, 0.0]).reshape(&[3, 3]);
let ortho = gram_schmidt(&matrix).expect("Gram-Schmidt should succeed");
let shape = ortho.shape();
let data = ortho.to_vec();
for j in 0..shape[1] {
let mut norm_sq = 0.0;
for i in 0..shape[0] {
norm_sq += data[i * shape[1] + j] * data[i * shape[1] + j];
}
assert!((norm_sq - 1.0).abs() < 1e-10, "Column norm should be 1");
}
}
#[test]
fn test_inverse_2x2() {
let matrix = Array::from_vec(vec![4.0, 7.0, 2.0, 6.0]).reshape(&[2, 2]);
let inv = inverse_2x2_or_general(&matrix).expect("matrix inverse should succeed");
let product = matrix_multiply(&matrix, &inv).expect("matrix multiply should succeed");
assert!((product.get(&[0, 0]).expect("valid index") - 1.0).abs() < 1e-10);
assert!((product.get(&[0, 1]).expect("valid index")).abs() < 1e-10);
assert!((product.get(&[1, 0]).expect("valid index")).abs() < 1e-10);
assert!((product.get(&[1, 1]).expect("valid index") - 1.0).abs() < 1e-10);
}
#[test]
fn test_nonnegative_cp_als() {
let tensor =
Array::from_vec(vec![1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0]).reshape(&[2, 2, 2]);
let result =
nonnegative_cp_als(&tensor, 2, 50, 1e-4).expect("nonnegative CP-ALS should succeed");
for factor in &result.factors {
for val in factor.to_vec() {
assert!(val >= 0.0, "All elements should be non-negative");
}
}
}
#[test]
fn test_tucker_reconstruct() {
let tensor =
Array::from_vec(vec![1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0]).reshape(&[2, 2, 2]);
let result =
tucker_decomposition(&tensor, &[2, 2, 2]).expect("tucker decomposition should succeed");
let reconstructed = tucker_reconstruct(&result).expect("tucker reconstruct should succeed");
assert_eq!(reconstructed.shape(), tensor.shape());
let error =
frobenius_error(&tensor, &reconstructed).expect("frobenius error should succeed");
let frobenius_norm: f64 = tensor.to_vec().iter().map(|x| x * x).sum::<f64>().sqrt();
let relative_error = error / frobenius_norm;
assert!(
relative_error < 2.0, "Tucker relative error {} indicates algorithm failure",
relative_error
);
}
}