1use crate::error::SvdLibError;
2use crate::utils::determine_chunk_size;
3use crate::{Diagnostics, SMat, SvdFloat, SvdRec};
4use nalgebra_sparse::na::{ComplexField, DMatrix, DVector, RealField};
5use ndarray::Array1;
6use nshare::IntoNdarray2;
7use rand::prelude::{Distribution, StdRng};
8use rand::SeedableRng;
9use rand_distr::Normal;
10use rayon::iter::ParallelIterator;
11use rayon::prelude::{IndexedParallelIterator, IntoParallelIterator};
12use std::ops::Mul;
13use std::time::Instant;
14
15#[derive(Debug, Clone, Copy, PartialEq)]
16pub enum PowerIterationNormalizer {
17 QR,
18 LU,
19 None,
20}
21
22const PARALLEL_THRESHOLD_ROWS: usize = 5000;
23const PARALLEL_THRESHOLD_COLS: usize = 1000;
24const PARALLEL_THRESHOLD_ELEMENTS: usize = 100_000;
25
26pub fn randomized_svd<T, M>(
27 m: &M,
28 target_rank: usize,
29 n_oversamples: usize,
30 n_power_iters: usize,
31 power_iteration_normalizer: PowerIterationNormalizer,
32 mean_center: bool,
33 seed: Option<u64>,
34 verbose: bool,
35) -> anyhow::Result<SvdRec<T>>
36where
37 T: SvdFloat + RealField,
38 M: SMat<T> + std::marker::Sync,
39 T: ComplexField,
40{
41 let start = Instant::now();
42 let m_rows = m.nrows();
43 let m_cols = m.ncols();
44
45 let rank = target_rank.min(m_rows.min(m_cols));
46 let l = rank + n_oversamples;
47
48 let column_means: Option<DVector<T>> = if mean_center {
49 if verbose {
50 println!("SVD | Randomized | Computing column means....");
51 }
52 Some(DVector::from(m.compute_column_means()))
53 } else {
54 None
55 };
56 if verbose && mean_center {
57 println!(
58 "SVD | Randomized | Computed column means, took: {:?} of total running time",
59 start.elapsed()
60 );
61 }
62
63 let omega = generate_random_matrix(m_cols, l, seed);
64
65 let mut y = DMatrix::<T>::zeros(m_rows, l);
66 if verbose {
67 println!("SVD | Randomized | Multiplying m with omega matrix....");
68 }
69 multiply_matrix_centered(m, &omega, &mut y, false, &column_means);
70 if verbose {
71 println!(
72 "SVD | Randomized | Multiplication done, took: {:?} of total running time",
73 start.elapsed()
74 );
75 }
76 if verbose {
77 println!("SVD | Randomized | Starting power iterations....");
78 }
79 if n_power_iters > 0 {
80 let mut z = DMatrix::<T>::zeros(m_cols, l);
81
82 for i in 0..n_power_iters {
83 if verbose {
84 println!(
85 "SVD | Randomized | Running power-iteration: {:?}, current time: {:?}",
86 i,
87 start.elapsed()
88 );
89 }
90 multiply_matrix_centered(m, &y, &mut z, true, &column_means);
91 if verbose {
92 println!(
93 "SVD | Randomized | Forward Multiplication {:?}",
94 start.elapsed()
95 );
96 }
97 match power_iteration_normalizer {
98 PowerIterationNormalizer::QR => {
99 let qr = z.qr();
100 z = qr.q();
101 }
102 PowerIterationNormalizer::LU => {
103 normalize_columns(&mut z);
104 }
105 PowerIterationNormalizer::None => {}
106 }
107 if verbose {
108 println!(
109 "SVD | Randomized | Power Iteration Normalization Forward-Step {:?}",
110 start.elapsed()
111 );
112 }
113
114 multiply_matrix_centered(m, &z, &mut y, false, &column_means);
115 if verbose {
116 println!(
117 "SVD | Randomized | Backward Multiplication {:?}",
118 start.elapsed()
119 );
120 }
121 match power_iteration_normalizer {
122 PowerIterationNormalizer::QR => {
123 let qr = y.qr();
124 y = qr.q();
125 }
126 PowerIterationNormalizer::LU => normalize_columns(&mut y),
127 PowerIterationNormalizer::None => {}
128 }
129 if verbose {
130 println!(
131 "SVD | Randomized | Power Iteration Normalization Backward-Step {:?}",
132 start.elapsed()
133 );
134 }
135 }
136 }
137 if verbose {
138 println!(
139 "SVD | Randomized | Running QR-Normalization after Power-Iterations {:?}",
140 start.elapsed()
141 );
142 }
143 let qr = y.qr();
144 let y = qr.q();
145 if verbose {
146 println!(
147 "SVD | Randomized | Finished QR-Normalization after Power-Iterations {:?}",
148 start.elapsed()
149 );
150 }
151
152 let mut b = DMatrix::<T>::zeros(y.ncols(), m_cols);
153 multiply_transposed_by_matrix_centered(&y, m, &mut b, &column_means);
154 if verbose {
155 println!(
156 "SVD | Randomized | Transposed Matrix Multiplication {:?}",
157 start.elapsed()
158 );
159 }
160 let svd = b.svd(true, true);
161 if verbose {
162 println!(
163 "SVD | Randomized | Running Singular Value Decomposition, took {:?}",
164 start.elapsed()
165 );
166 }
167 let u_b = svd
168 .u
169 .ok_or_else(|| SvdLibError::Las2Error("SVD U computation failed".to_string()))?;
170 let singular_values = svd.singular_values;
171 let vt = svd
172 .v_t
173 .ok_or_else(|| SvdLibError::Las2Error("SVD V_t computation failed".to_string()))?;
174
175 let u = y.mul(&u_b);
176 let actual_rank = target_rank.min(singular_values.len());
177
178 let u_subset = u.columns(0, actual_rank);
179 let s = convert_singular_values(
180 <DVector<T>>::from(singular_values.rows(0, actual_rank)),
181 actual_rank,
182 );
183 let vt_subset = vt.rows(0, actual_rank).into_owned();
184 let u = u_subset.into_owned().into_ndarray2();
185 let vt = vt_subset.into_ndarray2();
186 Ok(SvdRec {
187 d: actual_rank,
188 u,
189 s,
190 vt,
191 diagnostics: create_diagnostics(
192 m,
193 actual_rank,
194 target_rank,
195 n_power_iters,
196 seed.unwrap_or(0) as u32,
197 ),
198 })
199}
200
201fn convert_singular_values<T: SvdFloat + ComplexField>(
202 values: DVector<T::RealField>,
203 size: usize,
204) -> Array1<T> {
205 let mut array = Array1::zeros(size);
206
207 for i in 0..size {
208 array[i] = T::from_real(values[i].clone());
209 }
210
211 array
212}
213
214fn compute_column_means<T, M>(m: &M) -> Option<DVector<T>>
215where
216 T: SvdFloat + RealField + Send + Sync,
217 M: SMat<T> + Sync,
218{
219 let m_rows = m.nrows();
220 let m_cols = m.ncols();
221
222 let means: Vec<T> = (0..m_cols)
223 .into_par_iter()
224 .map(|j| {
225 let mut col_vec = vec![T::zero(); m_cols];
226 let mut result_vec = vec![T::zero(); m_rows];
227
228 col_vec[j] = T::one();
229
230 m.svd_opa(&col_vec, &mut result_vec, false);
231
232 let sum: T = result_vec.iter().copied().sum();
233 sum / T::from_f64(m_rows as f64).unwrap()
234 })
235 .collect();
236
237 Some(DVector::from_vec(means))
238}
239
240fn create_diagnostics<T, M: SMat<T>>(
241 a: &M,
242 d: usize,
243 target_rank: usize,
244 power_iterations: usize,
245 seed: u32,
246) -> Diagnostics<T>
247where
248 T: SvdFloat,
249{
250 Diagnostics {
251 non_zero: a.nnz(),
252 dimensions: target_rank,
253 iterations: power_iterations,
254 transposed: false,
255 lanczos_steps: 0, ritz_values_stabilized: d,
257 significant_values: d,
258 singular_values: d,
259 end_interval: [T::from(-1e-30).unwrap(), T::from(1e-30).unwrap()],
260 kappa: T::from(1e-6).unwrap(),
261 random_seed: seed,
262 }
263}
264
265fn normalize_columns<T: SvdFloat + RealField + Send + Sync>(matrix: &mut DMatrix<T>) {
266 let rows = matrix.nrows();
267 let cols = matrix.ncols();
268
269 if rows < PARALLEL_THRESHOLD_ROWS && cols < PARALLEL_THRESHOLD_COLS {
270 for j in 0..cols {
271 let mut norm = T::zero();
272
273 for i in 0..rows {
275 norm += ComplexField::powi(matrix[(i, j)], 2);
276 }
277 norm = ComplexField::sqrt(norm);
278
279 if norm > T::from_f64(1e-10).unwrap() {
280 let scale = T::one() / norm;
281 for i in 0..rows {
282 matrix[(i, j)] *= scale;
283 }
284 }
285 }
286 return;
287 }
288
289 let norms: Vec<T> = (0..cols)
290 .into_par_iter()
291 .map(|j| {
292 let mut norm = T::zero();
293 for i in 0..rows {
294 let val = unsafe { *matrix.get_unchecked((i, j)) };
295 norm += ComplexField::powi(val, 2);
296 }
297 ComplexField::sqrt(norm)
298 })
299 .collect();
300
301 let scales: Vec<(usize, T)> = norms
302 .into_iter()
303 .enumerate()
304 .filter_map(|(j, norm)| {
305 if norm > T::from_f64(1e-10).unwrap() {
306 Some((j, T::one() / norm))
307 } else {
308 None }
310 })
311 .collect();
312
313 scales.iter().for_each(|(j, scale)| {
314 for i in 0..rows {
315 let value = matrix.get_mut((i, *j)).unwrap();
316 *value = value.clone() * scale.clone();
317 }
318 });
319}
320
321fn generate_random_matrix<T: SvdFloat + RealField>(
326 rows: usize,
327 cols: usize,
328 seed: Option<u64>,
329) -> DMatrix<T> {
330 let mut rng = match seed {
331 Some(s) => StdRng::seed_from_u64(s),
332 None => StdRng::seed_from_u64(0),
333 };
334
335 let normal = Normal::new(0.0, 1.0).unwrap();
336 DMatrix::from_fn(rows, cols, |_, _| {
337 T::from_f64(normal.sample(&mut rng)).unwrap()
338 })
339}
340
341fn multiply_matrix<T: SvdFloat, M: SMat<T>>(
342 sparse: &M,
343 dense: &DMatrix<T>,
344 result: &mut DMatrix<T>,
345 transpose_sparse: bool,
346) {
347 sparse.multiply_with_dense(dense, result, transpose_sparse)
348}
349
350fn multiply_transposed_by_matrix<T: SvdFloat, M: SMat<T> + std::marker::Sync>(
351 q: &DMatrix<T>,
352 sparse: &M,
353 result: &mut DMatrix<T>,
354) {
355 let q_rows = q.nrows();
356 let q_cols = q.ncols();
357 let sparse_rows = sparse.nrows();
358 let sparse_cols = sparse.ncols();
359
360 eprintln!("Q dimensions: {} x {}", q_rows, q_cols);
361 eprintln!("Sparse dimensions: {} x {}", sparse_rows, sparse_cols);
362 eprintln!("Result dimensions: {} x {}", result.nrows(), result.ncols());
363
364 assert_eq!(
365 q_rows, sparse_rows,
366 "Dimension mismatch: Q has {} rows but sparse has {} rows",
367 q_rows, sparse_rows
368 );
369
370 assert_eq!(
371 result.nrows(),
372 q_cols,
373 "Result matrix has incorrect row count: expected {}, got {}",
374 q_cols,
375 result.nrows()
376 );
377 assert_eq!(
378 result.ncols(),
379 sparse_cols,
380 "Result matrix has incorrect column count: expected {}, got {}",
381 sparse_cols,
382 result.ncols()
383 );
384
385 let chunk_size = determine_chunk_size(q_cols);
386
387 let chunk_results: Vec<Vec<(usize, Vec<T>)>> = (0..q_cols)
388 .into_par_iter()
389 .chunks(chunk_size)
390 .map(|chunk| {
391 let mut chunk_results = Vec::with_capacity(chunk.len());
392
393 for &col_idx in &chunk {
394 let mut q_col = vec![T::zero(); q_rows];
395 for i in 0..q_rows {
396 q_col[i] = q[(i, col_idx)];
397 }
398
399 let mut result_row = vec![T::zero(); sparse_cols];
400
401 sparse.svd_opa(&q_col, &mut result_row, true);
402
403 chunk_results.push((col_idx, result_row));
404 }
405 chunk_results
406 })
407 .collect();
408
409 for chunk_result in chunk_results {
410 for (row_idx, row_values) in chunk_result {
411 for j in 0..sparse_cols {
412 result[(row_idx, j)] = row_values[j];
413 }
414 }
415 }
416}
417
418pub fn svd_flip<T: SvdFloat + 'static>(
419 u: Option<&mut DMatrix<T>>,
420 v: Option<&mut DMatrix<T>>,
421 u_based_decision: bool,
422) -> Result<(), SvdLibError> {
423 if u.is_none() && v.is_none() {
424 return Err(SvdLibError::Las2Error(
425 "Both u and v cannot be None".to_string(),
426 ));
427 }
428
429 if u_based_decision {
430 if u.is_none() {
431 return Err(SvdLibError::Las2Error(
432 "u cannot be None when u_based_decision is true".to_string(),
433 ));
434 }
435
436 let u = u.unwrap();
437 let ncols = u.ncols();
438 let nrows = u.nrows();
439
440 let mut signs = DVector::from_element(ncols, T::one());
441
442 for j in 0..ncols {
443 let mut max_abs = T::zero();
444 let mut max_idx = 0;
445
446 for i in 0..nrows {
447 let abs_val = num_traits::Float::abs(u[(i, j)]);
448 if abs_val > max_abs {
449 max_abs = abs_val;
450 max_idx = i;
451 }
452 }
453
454 if u[(max_idx, j)] < T::zero() {
455 signs[j] = -T::one();
456 }
457 }
458
459 for j in 0..ncols {
460 for i in 0..nrows {
461 u[(i, j)] *= signs[j];
462 }
463 }
464
465 if let Some(v) = v {
466 let v_nrows = v.nrows();
467 let v_ncols = v.ncols();
468
469 for i in 0..v_nrows.min(signs.len()) {
470 for j in 0..v_ncols {
471 v[(i, j)] *= signs[i];
472 }
473 }
474 }
475 } else {
476 if v.is_none() {
477 return Err(SvdLibError::Las2Error(
478 "v cannot be None when u_based_decision is false".to_string(),
479 ));
480 }
481
482 let v = v.unwrap();
483 let nrows = v.nrows();
484 let ncols = v.ncols();
485
486 let mut signs = DVector::from_element(nrows, T::one());
487
488 for i in 0..nrows {
489 let mut max_abs = T::zero();
490 let mut max_idx = 0;
491
492 for j in 0..ncols {
493 let abs_val = num_traits::Float::abs(v[(i, j)]);
494 if abs_val > max_abs {
495 max_abs = abs_val;
496 max_idx = j;
497 }
498 }
499
500 if v[(i, max_idx)] < T::zero() {
501 signs[i] = -T::one();
502 }
503 }
504
505 for i in 0..nrows {
506 for j in 0..ncols {
507 v[(i, j)] *= signs[i];
508 }
509 }
510
511 if let Some(u) = u {
512 let u_nrows = u.nrows();
513 let u_ncols = u.ncols();
514
515 for j in 0..u_ncols.min(signs.len()) {
516 for i in 0..u_nrows {
517 u[(i, j)] *= signs[j];
518 }
519 }
520 }
521 }
522
523 Ok(())
524}
525
526fn multiply_matrix_centered<T: SvdFloat, M: SMat<T> + std::marker::Sync>(
527 sparse: &M,
528 dense: &DMatrix<T>,
529 result: &mut DMatrix<T>,
530 transpose_sparse: bool,
531 column_means: &Option<DVector<T>>,
532) {
533 if column_means.is_none() {
534 multiply_matrix(sparse, dense, result, transpose_sparse);
535 return;
536 }
537
538 let means = column_means.as_ref().unwrap();
539 sparse.multiply_with_dense_centered(dense, result, transpose_sparse, means)
540}
541
542fn multiply_transposed_by_matrix_centered<T: SvdFloat, M: SMat<T> + std::marker::Sync>(
543 q: &DMatrix<T>,
544 sparse: &M,
545 result: &mut DMatrix<T>,
546 column_means: &Option<DVector<T>>,
547) {
548 if column_means.is_none() {
550 multiply_transposed_by_matrix(q, sparse, result);
551 return;
552 }
553
554 let means = column_means.as_ref().unwrap();
555 let q_rows = q.nrows();
556 let q_cols = q.ncols();
557 let sparse_rows = sparse.nrows();
558 let sparse_cols = sparse.ncols();
559
560 assert_eq!(
561 q_rows, sparse_rows,
562 "Dimension mismatch: Q has {} rows but sparse has {} rows",
563 q_rows, sparse_rows
564 );
565
566 assert_eq!(
567 result.nrows(),
568 q_cols,
569 "Result matrix has incorrect row count: expected {}, got {}",
570 q_cols,
571 result.nrows()
572 );
573 assert_eq!(
574 result.ncols(),
575 sparse_cols,
576 "Result matrix has incorrect column count: expected {}, got {}",
577 sparse_cols,
578 result.ncols()
579 );
580
581 let chunk_size = determine_chunk_size(q_cols);
582
583 let chunk_results: Vec<Vec<(usize, Vec<T>)>> = (0..q_cols)
584 .into_par_iter()
585 .chunks(chunk_size)
586 .map(|chunk| {
587 let mut chunk_results = Vec::with_capacity(chunk.len());
588
589 for &col_idx in &chunk {
590 let mut q_col = vec![T::zero(); q_rows];
591 for i in 0..q_rows {
592 q_col[i] = q[(i, col_idx)];
593 }
594
595 let mut result_row = vec![T::zero(); sparse_cols];
596
597 sparse.svd_opa(&q_col, &mut result_row, true);
598
599 let mut q_sum = T::zero();
600 for &val in &q_col {
601 q_sum += val;
602 }
603
604 for j in 0..sparse_cols {
605 result_row[j] -= means[j] * q_sum;
606 }
607
608 chunk_results.push((col_idx, result_row));
609 }
610 chunk_results
611 })
612 .collect();
613
614 for chunk_result in chunk_results {
615 for (row_idx, row_values) in chunk_result {
616 for j in 0..sparse_cols {
617 result[(row_idx, j)] = row_values[j];
618 }
619 }
620 }
621}
622
623#[cfg(test)]
624mod randomized_svd_tests {
625 use super::*;
626 use crate::randomized::{randomized_svd, PowerIterationNormalizer};
627 use nalgebra_sparse::coo::CooMatrix;
628 use nalgebra_sparse::CsrMatrix;
629 use ndarray::Array2;
630 use rand::rngs::StdRng;
631 use rand::{Rng, SeedableRng};
632 use rayon::ThreadPoolBuilder;
633 use std::sync::Once;
634
635 static INIT: Once = Once::new();
636
637 fn setup_thread_pool() {
638 INIT.call_once(|| {
639 ThreadPoolBuilder::new()
640 .num_threads(16)
641 .build_global()
642 .expect("Failed to build global thread pool");
643
644 println!("Initialized thread pool with {} threads", 16);
645 });
646 }
647
648 fn create_sparse_matrix(
649 rows: usize,
650 cols: usize,
651 density: f64,
652 ) -> nalgebra_sparse::coo::CooMatrix<f64> {
653 use std::collections::HashSet;
654
655 let mut coo = nalgebra_sparse::coo::CooMatrix::new(rows, cols);
656
657 let mut rng = StdRng::seed_from_u64(42);
658
659 let nnz = (rows as f64 * cols as f64 * density).round() as usize;
660
661 let nnz = nnz.max(1);
662
663 let mut positions = HashSet::new();
664
665 while positions.len() < nnz {
666 let i = rng.gen_range(0..rows);
667 let j = rng.gen_range(0..cols);
668
669 if positions.insert((i, j)) {
670 let val = loop {
671 let v: f64 = rng.gen_range(-10.0..10.0);
672 if v.abs() > 1e-10 {
673 break v;
674 }
675 };
676
677 coo.push(i, j, val);
678 }
679 }
680
681 let actual_density = coo.nnz() as f64 / (rows as f64 * cols as f64);
682 println!("Created sparse matrix: {} x {}", rows, cols);
683 println!(" - Requested density: {:.6}", density);
684 println!(" - Actual density: {:.6}", actual_density);
685 println!(" - Sparsity: {:.4}%", (1.0 - actual_density) * 100.0);
686 println!(" - Non-zeros: {}", coo.nnz());
687
688 coo
689 }
690
691 #[test]
692 fn test_randomized_svd_accuracy() {
693 setup_thread_pool();
694
695 let coo = create_sparse_matrix(500, 40, 0.1);
696
697
698 let csr = CsrMatrix::from(&coo);
699
700 let mut std_svd = crate::lanczos::svd_dim_seed(&csr, 10, 42).unwrap();
701
702 let rand_svd = randomized_svd(
703 &csr,
704 10,
705 5,
706 4,
707 PowerIterationNormalizer::QR,
708 false,
709 Some(42),
710 true,
711 )
712 .unwrap();
713
714 assert_eq!(rand_svd.d, 10, "Expected rank of 10");
715
716 let rel_tol = 0.4;
717 let compare_count = std::cmp::min(std_svd.d, rand_svd.d);
718 println!("Standard SVD has {} dimensions", std_svd.d);
719 println!("Randomized SVD has {} dimensions", rand_svd.d);
720
721 for i in 0..compare_count {
722 let rel_diff = (std_svd.s[i] - rand_svd.s[i]).abs() / std_svd.s[i];
723 println!(
724 "Singular value {}: standard={}, randomized={}, rel_diff={}",
725 i, std_svd.s[i], rand_svd.s[i], rel_diff
726 );
727 assert!(
728 rel_diff < rel_tol,
729 "Dominant singular value {} differs too much: rel diff = {}, standard = {}, randomized = {}",
730 i, rel_diff, std_svd.s[i], rand_svd.s[i]
731 );
732 }
733
734
735 }
736
737 #[test]
739 fn test_randomized_svd_with_mean_centering() {
740 setup_thread_pool();
741
742 let mut coo = CooMatrix::<f64>::new(30, 10);
743 let mut rng = StdRng::seed_from_u64(123);
744
745 let column_means: Vec<f64> = (0..10).map(|i| i as f64 * 2.0).collect();
746
747 let mut u = vec![vec![0.0; 3]; 30]; let mut v = vec![vec![0.0; 3]; 10];
749
750 for i in 0..30 {
751 for j in 0..3 {
752 u[i][j] = rng.gen_range(-1.0..1.0);
753 }
754 }
755
756 for i in 0..10 {
757 for j in 0..3 {
758 v[i][j] = rng.gen_range(-1.0..1.0);
759 }
760 }
761
762 for i in 0..30 {
763 for j in 0..10 {
764 let mut val = 0.0;
765 for k in 0..3 {
766 val += u[i][k] * v[j][k];
767 }
768 val = val + column_means[j] + rng.gen_range(-0.1..0.1);
769 coo.push(i, j, val);
770 }
771 }
772
773 let csr = CsrMatrix::from(&coo);
774
775 let svd_no_center = randomized_svd(
776 &csr,
777 3,
778 3,
779 2,
780 PowerIterationNormalizer::QR,
781 false,
782 Some(42),
783 false,
784 )
785 .unwrap();
786
787 let svd_with_center = randomized_svd(
788 &csr,
789 3,
790 3,
791 2,
792 PowerIterationNormalizer::QR,
793 true,
794 Some(42),
795 false,
796 )
797 .unwrap();
798
799 println!("Singular values without centering: {:?}", svd_no_center.s);
800 println!("Singular values with centering: {:?}", svd_with_center.s);
801 }
802
803 #[test]
804 fn test_randomized_svd_large_sparse() {
805 setup_thread_pool();
806
807 let test_matrix = create_sparse_matrix(5000, 1000, 0.01);
808
809 let csr = CsrMatrix::from(&test_matrix);
810
811 let result = randomized_svd(
812 &csr,
813 20,
814 10,
815 2,
816 PowerIterationNormalizer::QR,
817 false,
818 Some(42),
819 false,
820 );
821
822 assert!(
823 result.is_ok(),
824 "Randomized SVD failed on large sparse matrix: {:?}",
825 result.err().unwrap()
826 );
827
828 let svd = result.unwrap();
829 assert_eq!(svd.d, 20, "Expected rank of 20");
830 assert_eq!(svd.u.ncols(), 20, "Expected 20 left singular vectors");
831 assert_eq!(svd.u.nrows(), 5000, "Expected 5000 columns in U transpose");
832 assert_eq!(svd.vt.nrows(), 20, "Expected 20 right singular vectors");
833 assert_eq!(svd.vt.ncols(), 1000, "Expected 1000 columns in V transpose");
834
835 for i in 1..svd.s.len() {
836 assert!(svd.s[i] > 0.0, "Singular values should be positive");
837 assert!(
838 svd.s[i - 1] >= svd.s[i],
839 "Singular values should be in descending order"
840 );
841 }
842 }
843
844 #[test]
846 fn test_power_iteration_impact() {
847 setup_thread_pool();
848
849 let mut coo = CooMatrix::<f64>::new(100, 50);
850 let mut rng = StdRng::seed_from_u64(987);
851
852 let mut u = vec![vec![0.0; 10]; 100];
853 let mut v = vec![vec![0.0; 10]; 50];
854
855 for i in 0..100 {
856 for j in 0..10 {
857 u[i][j] = rng.gen_range(-1.0..1.0);
858 }
859 }
860
861 for i in 0..50 {
862 for j in 0..10 {
863 v[i][j] = rng.gen_range(-1.0..1.0);
864 }
865 }
866
867 for i in 0..100 {
868 for j in 0..50 {
869 let mut val = 0.0;
870 for k in 0..10 {
871 val += u[i][k] * v[j][k];
872 }
873 val += rng.gen_range(-0.01..0.01);
874 coo.push(i, j, val);
875 }
876 }
877
878 let csr = CsrMatrix::from(&coo);
879
880 let powers = [0, 1, 3, 5];
881 let mut errors = Vec::new();
882
883 let mut dense_mat = Array2::<f64>::zeros((100, 50));
884 for (i, j, val) in csr.triplet_iter() {
885 dense_mat[[i, j]] = *val;
886 }
887 let matrix_norm = dense_mat.iter().map(|x| x.powi(2)).sum::<f64>().sqrt();
888
889 for &power in &powers {
890 let svd = randomized_svd(
891 &csr,
892 10,
893 5,
894 power,
895 PowerIterationNormalizer::QR,
896 false,
897 Some(42),
898 false,
899 )
900 .unwrap();
901
902 let recon = svd.recompose();
903 let mut error = 0.0;
904
905 for i in 0..100 {
906 for j in 0..50 {
907 error += (dense_mat[[i, j]] - recon[[i, j]]).powi(2);
908 }
909 }
910
911 error = error.sqrt() / matrix_norm;
912 errors.push(error);
913
914 println!("Power iterations: {}, Relative error: {}", power, error);
915 }
916
917 let mut improved = false;
918 for i in 1..errors.len() {
919 if errors[i] < errors[0] * 0.9 {
920 improved = true;
921 break;
922 }
923 }
924
925 assert!(
926 improved,
927 "Power iterations did not improve accuracy as expected"
928 );
929 }
930}