use scirs2_core::ndarray::{s, Array1, Array2, ArrayView2};
use scirs2_core::numeric::{Float, NumAssign, One, Zero};
use scirs2_core::random::{self, Rng, RngExt};
use std::fmt::Debug;
use std::iter::Sum;
use scirs2_core::validation::{check_2d, check_positive};
use crate::decomposition::{qr, svd};
use crate::error::{LinalgError, LinalgResult};
#[allow(dead_code)]
pub fn nmf<F>(
a: &ArrayView2<F>,
rank: usize,
max_iter: usize,
tol: F,
) -> LinalgResult<(Array2<F>, Array2<F>)>
where
F: Float + NumAssign + Zero + One + Sum + Debug + 'static + std::fmt::Display,
{
check_2d(a, "a")?;
check_positive(F::from(rank).expect("Operation failed"), "rank")?;
let (m, n) = (a.nrows(), a.ncols());
for i in 0..m {
for j in 0..n {
if a[[i, j]] < F::zero() {
return Err(LinalgError::InvalidInputError(
"Input matrix must be non-negative for NMF".to_string(),
));
}
}
}
if rank > m.min(n) {
return Err(LinalgError::InvalidInputError(format!(
"Rank must be less than or equal to min(rows, cols) = {}",
m.min(n)
)));
}
let epsilon = F::from(1e-5).expect("Operation failed");
let mut w = Array2::<F>::zeros((m, rank));
let mut h = Array2::<F>::zeros((rank, n));
let mut rng = scirs2_core::random::rng();
for i in 0..m {
for j in 0..rank {
w[[i, j]] = F::from(rng.random::<f64>()).expect("Operation failed") + epsilon;
}
}
for i in 0..rank {
for j in 0..n {
h[[i, j]] = F::from(rng.random::<f64>()).expect("Operation failed") + epsilon;
}
}
let mut prev_error = F::infinity();
for _ in 0..max_iter {
let wt = w.t();
let wt_a = wt.dot(a);
let wt_w = wt.dot(&w);
let wt_w_h = wt_w.dot(&h);
for i in 0..rank {
for j in 0..n {
let numerator = wt_a[[i, j]];
let denominator = wt_w_h[[i, j]];
if denominator > epsilon {
h[[i, j]] = h[[i, j]] * numerator / denominator;
}
}
}
let ht = h.t();
let a_ht = a.dot(&ht);
let w_h = w.dot(&h);
let w_h_ht = w_h.dot(&ht);
for i in 0..m {
for j in 0..rank {
let numerator = a_ht[[i, j]];
let denominator = w_h_ht[[i, j]];
if denominator > epsilon {
w[[i, j]] = w[[i, j]] * numerator / denominator;
}
}
}
let a_approx = w.dot(&h);
let mut error = F::zero();
for i in 0..m {
for j in 0..n {
let diff = a[[i, j]] - a_approx[[i, j]];
error += diff * diff;
}
}
error = error.sqrt();
if (prev_error - error).abs() < tol {
break;
}
prev_error = error;
}
Ok((w, h))
}
#[allow(dead_code)]
pub fn interpolative_decomposition<F>(
a: &ArrayView2<F>,
k: usize,
method: &str,
) -> LinalgResult<(Array2<F>, Array2<F>)>
where
F: Float
+ NumAssign
+ Zero
+ One
+ Sum
+ Debug
+ 'static
+ scirs2_core::ndarray::ScalarOperand
+ Send
+ Sync,
{
check_2d(a, "a")?;
let (m, n) = (a.nrows(), a.ncols());
if k > n || k == 0 {
return Err(LinalgError::InvalidInputError(format!(
"k must be between 1 and n (number of columns) = {n}"
)));
}
match method.to_lowercase().as_str() {
"qr" => {
let mut a_copy = a.to_owned();
let mut col_indices = Vec::with_capacity(k);
for i in 0..k {
let mut max_norm = F::zero();
let mut max_col = i;
for j in i..n {
let col = a_copy.column(j);
let norm = col.iter().fold(F::zero(), |acc, &x| acc + x * x).sqrt();
if norm > max_norm {
max_norm = norm;
max_col = j;
}
}
if max_col != i {
for row in 0..m {
let temp = a_copy[[row, i]];
a_copy[[row, i]] = a_copy[[row, max_col]];
a_copy[[row, max_col]] = temp;
}
col_indices.push(max_col);
} else {
col_indices.push(i);
}
if i < k - 1 && i < m {
let pivot = a_copy.column(i).to_owned();
let pivot_norm = pivot.iter().fold(F::zero(), |acc, &x| acc + x * x).sqrt();
if pivot_norm > F::epsilon() {
for j in (i + 1)..n {
let col = a_copy.column(j).to_owned();
let dot_product = pivot
.iter()
.zip(col.iter())
.fold(F::zero(), |acc, (&p, &c)| acc + p * c)
/ pivot_norm;
for row in 0..m {
a_copy[[row, j]] =
a_copy[[row, j]] - dot_product * a_copy[[row, i]] / pivot_norm;
}
}
}
}
}
let mut c = Array2::<F>::zeros((m, k));
for (i, &col_idx) in col_indices.iter().enumerate() {
for row in 0..m {
c[[row, i]] = a[[row, col_idx]];
}
}
let mut z = Array2::<F>::zeros((k, n));
for (i, &col_idx) in col_indices.iter().enumerate() {
for j in 0..n {
if j == col_idx {
z[[i, j]] = F::one();
} else {
let a_j = a.column(j).to_owned();
let c_view = c.view();
let ct = c.t();
let ctc = ct.dot(&c_view);
let cta = ct.dot(&a_j.view());
let (u, s, vt) = svd(&ctc.view(), false, None)?;
let mut s_inv = s.clone();
for si in s_inv.iter_mut() {
if *si > F::epsilon() {
*si = F::one() / *si;
} else {
*si = F::zero();
}
}
let vtrans = vt.t();
let utrans = u.t();
let temp1 = utrans.dot(&cta.view());
let mut temp2 = Array1::<F>::zeros(k);
for j in 0..k {
temp2[j] = s_inv[j] * temp1[j];
}
let coeffs = vtrans.dot(&temp2.view());
for coef_idx in 0..k {
z[[coef_idx, j]] = coeffs[coef_idx];
}
}
}
}
Ok((c, z))
}
"svd" => {
let (u, s, vt) = svd(a, false, None)?;
let u_k = u.slice(s![.., ..k]).to_owned();
let s_k = s.slice(s![..k]).to_owned();
let vt_k = vt.slice(s![..k, ..]).to_owned();
let mut column_scores = vec![F::zero(); n];
let v = vt_k.t();
for j in 0..n {
for i in 0..k {
column_scores[j] += v[[j, i]].powi(2) * s_k[i];
}
}
let mut indices: Vec<usize> = (0..n).collect();
indices.sort_by(|&a, &b| {
column_scores[b]
.partial_cmp(&column_scores[a])
.unwrap_or(std::cmp::Ordering::Equal)
});
let col_indices: Vec<usize> = indices.into_iter().take(k).collect();
let mut c = Array2::<F>::zeros((m, k));
for (i, &col_idx) in col_indices.iter().enumerate() {
for row in 0..m {
c[[row, i]] = a[[row, col_idx]];
}
}
let c_pinv = c.t().dot(&u_k);
let mut s_inv_diag = Array2::<F>::zeros((k, k));
for i in 0..k {
if s_k[i] > F::epsilon() {
s_inv_diag[[i, i]] = F::one() / s_k[i];
}
}
let temp = c_pinv.dot(&s_inv_diag);
let z = temp.dot(&vt_k);
Ok((c, z))
}
_ => Err(LinalgError::InvalidInputError(format!(
"Unknown method: {method}. Expected 'qr' or 'svd'"
))),
}
}
#[allow(dead_code)]
pub fn cur_decomposition<F>(
a: &ArrayView2<F>,
k: usize,
c_samples: Option<usize>,
r_samples: Option<usize>,
method: &str,
) -> LinalgResult<(Array2<F>, Array2<F>, Array2<F>)>
where
F: Float
+ NumAssign
+ Zero
+ One
+ Sum
+ Debug
+ 'static
+ scirs2_core::ndarray::ScalarOperand
+ Send
+ Sync,
{
check_2d(a, "a")?;
let (m, n) = (a.nrows(), a.ncols());
if k > m.min(n) || k == 0 {
return Err(LinalgError::InvalidInputError(format!(
"k must be between 1 and min(rows, cols) = {}",
m.min(n)
)));
}
let c_samples = c_samples.unwrap_or(2 * k);
let r_samples = r_samples.unwrap_or(2 * k);
if c_samples > n || r_samples > m {
return Err(LinalgError::InvalidInputError(
"Number of _samples cannot exceed matrix dimensions".to_string(),
));
}
match method.to_lowercase().as_str() {
"uniform" => {
let mut col_indices = Vec::with_capacity(c_samples);
let mut row_indices = Vec::with_capacity(r_samples);
while col_indices.len() < c_samples {
let idx = scirs2_core::random::rng().random_range(0..n);
if !col_indices.contains(&idx) {
col_indices.push(idx);
}
}
while row_indices.len() < r_samples {
let idx = scirs2_core::random::rng().random_range(0..m);
if !row_indices.contains(&idx) {
row_indices.push(idx);
}
}
let mut c = Array2::<F>::zeros((m, c_samples));
let mut r = Array2::<F>::zeros((r_samples, n));
for (c_idx, &col) in col_indices.iter().enumerate() {
for i in 0..m {
c[[i, c_idx]] = a[[i, col]];
}
}
for (r_idx, &row) in row_indices.iter().enumerate() {
for j in 0..n {
r[[r_idx, j]] = a[[row, j]];
}
}
let mut w = Array2::<F>::zeros((r_samples, c_samples));
for (r_idx, &row) in row_indices.iter().enumerate() {
for (c_idx, &col) in col_indices.iter().enumerate() {
w[[r_idx, c_idx]] = a[[row, col]];
}
}
let (u_w, s_w, vt_w) = svd(&w.view(), true, None)?;
let mut effective_rank = 0;
for i in 0..s_w.len() {
if s_w[i] > F::epsilon() * s_w[0] {
effective_rank += 1;
} else {
break;
}
}
let u_w_k = u_w.slice(s![.., ..effective_rank]).to_owned();
let vt_w_k = vt_w.slice(s![..effective_rank, ..]).to_owned();
let mut s_w_inv = Array2::<F>::zeros((effective_rank, effective_rank));
for i in 0..effective_rank {
if s_w[i] > F::epsilon() {
s_w_inv[[i, i]] = F::one() / s_w[i];
}
}
let v_w_k = vt_w_k.t();
let u_w_k_t = u_w_k.t();
let temp = v_w_k.dot(&s_w_inv);
let u = temp.dot(&u_w_k_t);
Ok((c, u, r))
}
"leverage" => {
let mut rng = scirs2_core::random::rng();
let omega = Array2::<F>::from_shape_fn((n, k + 5), |_| {
F::from(rng.random::<f64>() * 2.0 - 1.0).expect("Operation failed")
});
let y = a.dot(&omega);
let (q, _) = qr(&y.view(), None)?;
let qt = q.t();
let b = qt.dot(a);
let (_, s, vt) = svd(&b.view(), false, None)?;
let s_k = s.slice(s![..k]).to_owned();
let vt_k = vt.slice(s![..k, ..]).to_owned();
let mut col_leverage = Array1::<F>::zeros(n);
for j in 0..n {
for i in 0..k {
col_leverage[j] += vt_k[[i, j]] * vt_k[[i, j]];
}
}
let v_k = vt_k.t();
let mut s_inv = Array2::<F>::zeros((k, k));
for i in 0..k {
if s_k[i] > F::epsilon() {
s_inv[[i, i]] = F::one() / s_k[i];
}
}
let v_s_inv = v_k.dot(&s_inv);
let u_approx = a.dot(&v_s_inv);
let mut row_leverage = Array1::<F>::zeros(m);
for i in 0..m {
for j in 0..k {
row_leverage[i] += u_approx[[i, j]] * u_approx[[i, j]];
}
}
let mut col_indices = Vec::with_capacity(c_samples);
let mut row_indices = Vec::with_capacity(r_samples);
let col_sum = col_leverage.sum();
let row_sum = row_leverage.sum();
for j in 0..n {
col_leverage[j] /= col_sum;
}
for i in 0..m {
row_leverage[i] /= row_sum;
}
for _ in 0..c_samples {
let rand_val = F::from(rng.random::<f64>()).expect("Operation failed");
let mut cumsum = F::zero();
let mut selected = 0;
for (j, &prob) in col_leverage.iter().enumerate() {
cumsum += prob;
if rand_val <= cumsum {
selected = j;
break;
}
}
col_indices.push(selected);
}
for _ in 0..r_samples {
let rand_val = F::from(rng.random::<f64>()).expect("Operation failed");
let mut cumsum = F::zero();
let mut selected = 0;
for (i, &prob) in row_leverage.iter().enumerate() {
cumsum += prob;
if rand_val <= cumsum {
selected = i;
break;
}
}
row_indices.push(selected);
}
let mut c = Array2::<F>::zeros((m, c_samples));
let mut r = Array2::<F>::zeros((r_samples, n));
for (c_idx, &col) in col_indices.iter().enumerate() {
let scale = F::one()
/ (F::from(c_samples).expect("Operation failed") * col_leverage[col]).sqrt();
for i in 0..m {
c[[i, c_idx]] = a[[i, col]] * scale;
}
}
for (r_idx, &row) in row_indices.iter().enumerate() {
let scale = F::one()
/ (F::from(r_samples).expect("Operation failed") * row_leverage[row]).sqrt();
for j in 0..n {
r[[r_idx, j]] = a[[row, j]] * scale;
}
}
let (c_u, c_s, c_vt) = svd(&c.view(), false, None)?;
let (r_u, r_s, r_vt) = svd(&r.view(), false, None)?;
let c_u_k = c_u.slice(s![.., ..k]).to_owned();
let c_vt_k = c_vt.slice(s![..k, ..]).to_owned();
let r_u_k = r_u.slice(s![.., ..k]).to_owned();
let r_vt_k = r_vt.slice(s![..k, ..]).to_owned();
let mut c_s_inv = Array2::<F>::zeros((k, k));
let mut r_s_inv = Array2::<F>::zeros((k, k));
for i in 0..k {
if c_s[i] > F::epsilon() {
c_s_inv[[i, i]] = F::one() / c_s[i];
}
if r_s[i] > F::epsilon() {
r_s_inv[[i, i]] = F::one() / r_s[i];
}
}
let c_v_k = c_vt_k.t();
let c_ut_k = c_u_k.t();
let c_pseudo = c_v_k.dot(&c_s_inv).dot(&c_ut_k);
let r_v_k = r_vt_k.t();
let r_ut_k = r_u_k.t();
let r_pseudo = r_v_k.dot(&r_s_inv).dot(&r_ut_k);
let temp = c_pseudo.dot(a);
let u = temp.dot(&r_pseudo);
Ok((c, u, r))
}
_ => Err(LinalgError::InvalidInputError(format!(
"Unknown method: {method}. Expected 'uniform' or 'leverage'"
))),
}
}
#[allow(dead_code)]
pub fn rank_revealing_qr<F>(
a: &ArrayView2<F>,
tol: F,
) -> LinalgResult<(Array2<F>, Array2<F>, Array2<F>)>
where
F: Float
+ NumAssign
+ Zero
+ One
+ Sum
+ Debug
+ 'static
+ scirs2_core::ndarray::ScalarOperand
+ Send
+ Sync,
{
check_2d(a, "a")?;
let (m, n) = (a.nrows(), a.ncols());
let min_dim = m.min(n);
let mut q = Array2::<F>::eye(m);
let mut r = a.to_owned();
let mut p = Array2::<F>::eye(n);
let mut col_norms = Vec::with_capacity(n);
for j in 0..n {
let col = r.column(j);
let norm_sq = col.iter().fold(F::zero(), |acc, &x| acc + x * x);
col_norms.push(norm_sq);
}
for k in 0..min_dim {
let mut max_norm = F::zero();
let mut max_col = k;
for (j, &norm) in col_norms.iter().enumerate().skip(k).take(n - k) {
if norm > max_norm {
max_norm = norm;
max_col = j;
}
}
if max_norm.sqrt() <= tol {
break;
}
if max_col != k {
for i in 0..m {
let temp = r[[i, k]];
r[[i, k]] = r[[i, max_col]];
r[[i, max_col]] = temp;
}
for i in 0..n {
let temp = p[[i, k]];
p[[i, k]] = p[[i, max_col]];
p[[i, max_col]] = temp;
}
col_norms.swap(k, max_col);
}
let mut x = Array1::<F>::zeros(m - k);
for i in k..m {
x[i - k] = r[[i, k]];
}
let x_norm = x.iter().fold(F::zero(), |acc, &val| acc + val * val).sqrt();
if x_norm > F::epsilon() {
let alpha = if x[0] >= F::zero() { -x_norm } else { x_norm };
let mut v = x.clone();
v[0] -= alpha;
let v_norm = v.iter().fold(F::zero(), |acc, &val| acc + val * val).sqrt();
if v_norm > F::epsilon() {
for i in 0..v.len() {
v[i] /= v_norm;
}
for j in k..n {
let mut r_col = Array1::<F>::zeros(m - k);
for i in k..m {
r_col[i - k] = r[[i, j]];
}
let dot_product = v
.iter()
.zip(r_col.iter())
.fold(F::zero(), |acc, (&vi, &ri)| acc + vi * ri);
for i in k..m {
r[[i, j]] -=
F::from(2.0).expect("Operation failed") * v[i - k] * dot_product;
}
}
for i in 0..m {
let mut q_row = Array1::<F>::zeros(m - k);
for j in k..m {
q_row[j - k] = q[[i, j]];
}
let dot_product = q_row
.iter()
.zip(v.iter())
.fold(F::zero(), |acc, (&qi, &vi)| acc + qi * vi);
for j in k..m {
q[[i, j]] -=
F::from(2.0).expect("Operation failed") * dot_product * v[j - k];
}
}
}
}
for j in k + 1..n {
col_norms[j] = F::zero();
for i in k + 1..m {
col_norms[j] += r[[i, j]] * r[[i, j]];
}
}
}
Ok((q, r, p))
}
#[allow(dead_code)]
pub fn utv_decomposition<F>(
a: &ArrayView2<F>,
variant: &str,
tol: F,
) -> LinalgResult<(Array2<F>, Array2<F>, Array2<F>)>
where
F: Float
+ NumAssign
+ Zero
+ One
+ Sum
+ Debug
+ 'static
+ scirs2_core::ndarray::ScalarOperand
+ Send
+ Sync,
{
check_2d(a, "a")?;
match variant.to_lowercase().as_str() {
"urv" => {
let (q, r, p) = rank_revealing_qr(a, tol)?;
Ok((q, r, p))
}
"utv" => {
let (m, n) = (a.nrows(), a.ncols());
let (q, r) = qr(a, None)?;
let mut rank = 0;
for i in 0..m.min(n) {
if r[[i, i]].abs() > tol {
rank += 1;
} else {
break;
}
}
if rank == 0 {
return Ok((Array2::eye(m), Array2::zeros((m, n)), Array2::eye(n)));
}
let r11 = r.slice(s![..rank, ..rank]).to_owned();
let (u11, s11, v11t) = svd(&r11.view(), true, None)?;
let mut s11_diag = Array2::zeros((rank, rank));
for i in 0..rank {
s11_diag[[i, i]] = s11[i];
}
let mut u_mid = Array2::zeros((m, m));
for i in 0..rank {
for j in 0..rank {
u_mid[[i, j]] = u11[[i, j]];
}
}
for i in rank..m {
u_mid[[i, i]] = F::one();
}
let v11 = v11t.t();
let mut v_mid = Array2::zeros((n, n));
for i in 0..rank {
for j in 0..rank {
v_mid[[i, j]] = v11[[i, j]];
}
}
for i in rank..n {
v_mid[[i, i]] = F::one();
}
let u = q.dot(&u_mid);
let mut t = Array2::zeros((m, n));
for i in 0..rank {
for j in 0..rank {
t[[i, j]] = s11_diag[[i, j]];
}
}
for i in 0..rank {
for j in rank..n {
let r_val = r[[i, j]];
let mut transformed = F::zero();
for k in 0..rank {
transformed += u11[[i, k]] * r_val;
}
t[[i, j]] = transformed;
}
}
let v = v_mid;
Ok((u, t, v))
}
_ => Err(LinalgError::InvalidInputError(format!(
"Unknown variant: {variant}. Expected 'urv' or 'utv'"
))),
}
}
#[cfg(test)]
mod tests {
use super::*;
use approx::assert_relative_eq;
use scirs2_core::ndarray::array;
#[test]
fn test_nmf_simple() {
let a = array![[1.0, 2.0, 3.0], [4.0, 5.0, 6.0], [7.0, 8.0, 9.0]];
let (w, h) = nmf(&a.view(), 2, 100, 1e-4).expect("Operation failed");
assert_eq!(w.shape(), &[3, 2]);
assert_eq!(h.shape(), &[2, 3]);
for i in 0..w.shape()[0] {
for j in 0..w.shape()[1] {
assert!(w[[i, j]] >= 0.0);
}
}
for i in 0..h.shape()[0] {
for j in 0..h.shape()[1] {
assert!(h[[i, j]] >= 0.0);
}
}
let wh = w.dot(&h);
let mut error = 0.0;
for i in 0..3 {
for j in 0..3 {
error += (a[[i, j]] - wh[[i, j]]).powi(2);
}
}
error = error.sqrt();
assert!(error / 9.0 < 1.0);
}
#[test]
fn test_interpolative_decomposition() {
let a = array![
[1.0, 2.0, 3.0, 4.0],
[5.0, 6.0, 7.0, 8.0],
[9.0, 10.0, 11.0, 12.0]
];
let method = "svd";
match interpolative_decomposition(&a.view(), 2, method) {
Ok((c, z)) => {
assert_eq!(c.nrows(), a.nrows());
assert!(z.nrows() <= a.ncols());
}
Err(_) => {
}
}
}
#[test]
fn test_cur_decomposition() {
let a = array![
[1.0, 2.0, 3.0, 4.0],
[5.0, 6.0, 7.0, 8.0],
[9.0, 10.0, 11.0, 12.0]
];
let (c, u, r) =
cur_decomposition(&a.view(), 2, Some(2), Some(2), "uniform").expect("Operation failed");
assert_eq!(c.shape(), &[3, 2]);
assert_eq!(u.shape(), &[2, 2]);
assert_eq!(r.shape(), &[2, 4]);
let _approx = c.dot(&u).dot(&r);
}
#[test]
fn test_rank_revealing_qr() {
let a = array![[1.0, 2.0, 3.0], [4.0, 5.0, 6.0], [7.0, 8.0, 9.0]];
let (q, r, p) = rank_revealing_qr(&a.view(), 1e-10).expect("Operation failed");
assert_eq!(q.shape(), &[3, 3]);
assert_eq!(r.shape(), &[3, 3]);
assert_eq!(p.shape(), &[3, 3]);
let qt = q.t();
let qtq = qt.dot(&q);
for i in 0..3 {
for j in 0..3 {
if i == j {
assert_relative_eq!(qtq[[i, j]], 1.0, epsilon = 1e-6);
} else {
assert_relative_eq!(qtq[[i, j]], 0.0, epsilon = 1e-6);
}
}
}
for i in 0..3 {
let row_sum: f64 = p.row(i).iter().map(|&x| x.abs()).sum();
let col_sum: f64 = p.column(i).iter().map(|&x| x.abs()).sum();
assert_relative_eq!(row_sum, 1.0, epsilon = 1e-6);
assert_relative_eq!(col_sum, 1.0, epsilon = 1e-6);
}
let qr = q.dot(&r);
let qrpt = qr.dot(&p.t());
for i in 0..3 {
for j in 0..3 {
assert_relative_eq!(qrpt[[i, j]], a[[i, j]], epsilon = 1e-6);
}
}
assert!(r[[0, 0]].abs() > 1e-6);
assert!(r[[1, 1]].abs() > 1e-6);
assert!(r[[2, 2]].abs() < 1e-6 || r[[2, 2]].abs() / r[[0, 0]].abs() < 1e-6);
}
#[test]
fn test_utv_decomposition() {
let a = array![[1.0, 2.0, 3.0], [4.0, 5.0, 6.0], [7.0, 8.0, 9.0]];
let (u, t, v) = utv_decomposition(&a.view(), "urv", 1e-10).expect("Operation failed");
assert_eq!(u.shape(), &[3, 3]);
assert_eq!(t.shape(), &[3, 3]);
assert_eq!(v.shape(), &[3, 3]);
let ut = u.t();
let utu = ut.dot(&u);
for i in 0..3 {
for j in 0..3 {
if i == j {
assert_relative_eq!(utu[[i, j]], 1.0, epsilon = 1e-6);
} else {
assert_relative_eq!(utu[[i, j]], 0.0, epsilon = 1e-6);
}
}
}
let vt = v.t();
let vtv = vt.dot(&v);
for i in 0..3 {
for j in 0..3 {
if i == j {
assert_relative_eq!(vtv[[i, j]], 1.0, epsilon = 1e-6);
} else {
assert_relative_eq!(vtv[[i, j]], 0.0, epsilon = 1e-6);
}
}
}
let ut_prod = u.dot(&t);
let utv = ut_prod.dot(&vt);
for i in 0..3 {
for j in 0..3 {
assert_relative_eq!(utv[[i, j]], a[[i, j]], epsilon = 1e-6);
}
}
}
}