Getting Started
This guide will help you quickly get up and running with single-svdlib, a high-performance Rust library for computing Singular Value Decomposition (SVD) on sparse matrices.
Installation
Add single-svdlib to your Rust project by including it in your Cargo.toml
:
[dependencies]
single-svdlib = "1.0.4"
The library also requires nalgebra-sparse
for sparse matrix representations:
[dependencies]
single-svdlib = "1.0.4"
nalgebra-sparse = "0.10.0"
Basic Usage
Let's walk through a simple example of using single-svdlib to compute the SVD of a sparse matrix.
Importing Required Dependencies
use nalgebra_sparse::{coo::CooMatrix, csr::CsrMatrix};
use single_svdlib::lanczos;
use ndarray::prelude::*;
Creating a Sparse Matrix
You can create a sparse matrix using the Coordinate (COO) format:
// Create a 3x3 matrix
let mut coo = CooMatrix::<f64>::new(3, 3);
// Add non-zero elements using (row, column, value)
coo.push(0, 0, 1.0); coo.push(0, 1, 16.0); coo.push(0, 2, 49.0);
coo.push(1, 0, 4.0); coo.push(1, 1, 25.0); coo.push(1, 2, 64.0);
coo.push(2, 0, 9.0); coo.push(2, 1, 36.0); coo.push(2, 2, 81.0);
// Convert to CSR format for better performance
let csr = CsrMatrix::from(&coo);
Computing SVD
single-svdlib offers several convenient methods to compute SVD:
// Compute SVD with all dimensions
let result = lanczos::svd(&csr);
// Or with a specific number of dimensions
let result = lanczos::svd_dim(&csr, 2);
// Or with a fixed random seed for reproducibility
let result = lanczos::svd_dim_seed(&csr, 2, 42);
// Unwrap the result
let svd = result.unwrap();
Using the SVD Results
Once you have computed the SVD, you can access the decomposition components:
// Access the singular values
let singular_values = &svd.s;
println!("Singular values: {:?}", singular_values);
// Access the left singular vectors (U)
let left_vectors = &svd.u;
// Access the right singular vectors transposed (V^T)
let right_vectors_transposed = &svd.vt;
// Get dimensionality (rank) of the decomposition
let rank = svd.d;
println!("Rank: {}", rank);
// Reconstruct the original matrix from the SVD components
let reconstructed = svd.recompose();
Handling Errors
The SVD functions return a Result
that you should properly handle:
match lanczos::svd(&csr) {
Ok(svd) => {
println!("SVD computed successfully!");
println!("Singular values: {:?}", svd.s);
},
Err(err) => {
eprintln!("Failed to compute SVD: {}", err);
}
}
Advanced Usage
Using Randomized SVD for Large Matrices
For very large sparse matrices, the randomized SVD algorithm often provides better performance:
use single_svdlib::randomized::{randomized_svd, PowerIterationNormalizer};
let svd = randomized_svd(
&csr,
10, // Target rank
5, // Oversampling parameter
3, // Number of power iterations
PowerIterationNormalizer::QR, // Normalization method
false, // Mean centering
Some(42), // Random seed
false // Verbose output
).unwrap();
Column Masking for Subspace SVD
If you only need to compute SVD on a subset of columns:
use single_svdlib::lanczos::masked::MaskedCSRMatrix;
// Specify which columns to include
let columns = vec![0, 2]; // Only use columns 0 and 2
// Create a masked view of the matrix
let masked_matrix = MaskedCSRMatrix::with_columns(&csr, &columns);
// Compute SVD on the masked matrix
let svd = lanczos::svd(&masked_matrix).unwrap();
Accessing Diagnostics
Examine computation details to understand algorithm behavior:
// Print diagnostic information
println!("Non-zero elements: {}", svd.diagnostics.non_zero);
println!("Algorithm iterations: {}", svd.diagnostics.iterations);
println!("Transposed during computation: {}", svd.diagnostics.transposed);
println!("Lanczos steps: {}", svd.diagnostics.lanczos_steps);
println!("Significant values found: {}", svd.diagnostics.significant_values);
println!("Random seed used: {}", svd.diagnostics.random_seed);
Complete Example
Here's a complete working example:
use nalgebra_sparse::{coo::CooMatrix, csr::CsrMatrix};
use single_svdlib::lanczos;
use std::error::Error;
fn main() -> Result<(), Box<dyn Error>> {
// Create a sparse matrix in COO format
let mut coo = CooMatrix::<f64>::new(3, 3);
coo.push(0, 0, 1.0); coo.push(0, 1, 16.0); coo.push(0, 2, 49.0);
coo.push(1, 0, 4.0); coo.push(1, 1, 25.0); coo.push(1, 2, 64.0);
coo.push(2, 0, 9.0); coo.push(2, 1, 36.0); coo.push(2, 2, 81.0);
// Convert to CSR for better performance
let csr = CsrMatrix::from(&coo);
// Compute SVD
let svd = lanczos::svd_dim_seed(&csr, 3, 42)?;
// Print the results
println!("Singular values:");
for (i, val) in svd.s.iter().enumerate() {
println!(" σ_{} = {:.6}", i, val);
}
// Print dimensions of the result matrices
println!("\nDimensions:");
println!(" U: {} × {}", svd.u.nrows(), svd.u.ncols());
println!(" Σ: length {}", svd.s.len());
println!(" V^T: {} × {}", svd.vt.nrows(), svd.vt.ncols());
// Reconstruct and verify
let reconstructed = svd.recompose();
println!("\nReconstructed matrix:");
for i in 0..3 {
for j in 0..3 {
print!("{:.1} ", reconstructed[[i, j]]);
}
println!();
}
Ok(())
}
Best Practices
Algorithm Selection
Choose the appropriate algorithm based on your matrix characteristics:
- Lanczos algorithm (
lanczos::svd_*
): Best for matrices up to ~10,000 × 10,000 with moderate sparsity - Randomized SVD (
randomized::randomized_svd
): Preferred for larger matrices or those with very high sparsity (>99%)
Performance Optimization
-
Matrix format matters:
- Convert COO matrices to CSR or CSC before computation
- CSR typically performs better for row-oriented operations
-
Tune randomized SVD parameters:
- Increase
n_power_iterations
for better accuracy (at the cost of speed) - Use
PowerIterationNormalizer::QR
for numerical stability - Adjust
n_oversamples
based on desired accuracy
- Increase
-
Memory efficiency:
- Request only the number of singular values you need using
svd_dim
- Use
f32
instead off64
when lower precision is acceptable
- Request only the number of singular values you need using
-
Parallel processing:
- The library uses Rayon for parallel processing by default
- For large matrices, ensure your system has adequate cores
Troubleshooting
Common Errors
-
Convergence issues: For extremely sparse matrices, try:
- Increasing the number of power iterations
- Using a larger value for
kappa
insvdLAS2
- Using randomized SVD instead of Lanczos
-
Dimension mismatch: Ensure that when using custom algorithm parameters that the requested dimensions don't exceed the matrix dimensions.
-
Memory limitations: For very large matrices:
- Compute only the dimensions you need
- Use column masking to focus on relevant subspaces
Debugging
Enable verbose output in randomized SVD to get more information:
let svd = randomized_svd(
&csr,
target_rank,
n_oversamples,
n_power_iterations,
PowerIterationNormalizer::QR,
false,
Some(42),
true // Set verbose to true
)?;
Next Steps
Now that you're familiar with the basics of single-svdlib, you might want to:
- Explore the API documentation for more detailed information
- Check out the examples directory in the repository
- Learn more about SVD applications in data analysis