pub struct Matrix { /* private fields */ }Expand description
A dense matrix stored in row-major order.
§Storage
Elements are stored contiguously: data[i * cols + j] holds A[i, j].
Implementations§
Source§impl Matrix
impl Matrix
Sourcepub fn identity(n: usize) -> Self
pub fn identity(n: usize) -> Self
Creates an identity matrix.
§Examples
use u_numflow::matrix::Matrix;
let eye = Matrix::identity(3);
assert_eq!(eye.get(0, 0), 1.0);
assert_eq!(eye.get(0, 1), 0.0);
assert_eq!(eye.get(2, 2), 1.0);Sourcepub fn transpose(&self) -> Self
pub fn transpose(&self) -> Self
Transpose: returns Aᵀ.
§Examples
use u_numflow::matrix::Matrix;
let m = Matrix::from_rows(&[&[1.0, 2.0, 3.0], &[4.0, 5.0, 6.0]]);
let t = m.transpose();
assert_eq!(t.rows(), 3);
assert_eq!(t.cols(), 2);
assert_eq!(t.get(0, 1), 4.0);Sourcepub fn add(&self, other: &Self) -> Result<Self, MatrixError>
pub fn add(&self, other: &Self) -> Result<Self, MatrixError>
Sourcepub fn sub(&self, other: &Self) -> Result<Self, MatrixError>
pub fn sub(&self, other: &Self) -> Result<Self, MatrixError>
Sourcepub fn mul_mat(&self, other: &Self) -> Result<Self, MatrixError>
pub fn mul_mat(&self, other: &Self) -> Result<Self, MatrixError>
Matrix multiplication: A · B.
Uses i-k-j loop order for better cache locality on row-major storage.
§Errors
Returns Err if self.cols != other.rows.
§Complexity
O(n·m·p) where self is n×m and other is m×p.
§Examples
use u_numflow::matrix::Matrix;
let a = Matrix::identity(3);
let b = Matrix::from_rows(&[&[1.0, 2.0, 3.0], &[4.0, 5.0, 6.0], &[7.0, 8.0, 9.0]]);
let c = a.mul_mat(&b).unwrap();
assert_eq!(c.get(2, 2), 9.0);Sourcepub fn frobenius_norm(&self) -> f64
pub fn frobenius_norm(&self) -> f64
Frobenius norm: ‖A‖_F = √(Σᵢⱼ aᵢⱼ²).
Sourcepub fn is_symmetric(&self, tol: f64) -> bool
pub fn is_symmetric(&self, tol: f64) -> bool
Checks whether the matrix is symmetric within tolerance.
Sourcepub fn determinant(&self) -> Result<f64, MatrixError>
pub fn determinant(&self) -> Result<f64, MatrixError>
Determinant via LU decomposition with partial pivoting.
§Errors
Returns Err(NotSquare) if the matrix is not square.
Returns Err(Singular) if a zero pivot is encountered.
§Complexity
O(n³/3).
§Examples
use u_numflow::matrix::Matrix;
let m = Matrix::from_rows(&[&[2.0, 3.0], &[1.0, 4.0]]);
assert!((m.determinant().unwrap() - 5.0).abs() < 1e-10);Sourcepub fn inverse(&self) -> Result<Self, MatrixError>
pub fn inverse(&self) -> Result<Self, MatrixError>
Matrix inverse via Gauss-Jordan elimination with partial pivoting.
§Algorithm
Augments [A | I], reduces to [I | A⁻¹] using row operations.
Reference: Golub & Van Loan (1996), Matrix Computations, §1.2.
§Errors
Returns Err(Singular) if the matrix is singular.
§Examples
use u_numflow::matrix::Matrix;
let a = Matrix::from_rows(&[&[4.0, 7.0], &[2.0, 6.0]]);
let inv = a.inverse().unwrap();
let eye = a.mul_mat(&inv).unwrap();
assert!((eye.get(0, 0) - 1.0).abs() < 1e-10);
assert!(eye.get(0, 1).abs() < 1e-10);Sourcepub fn cholesky(&self) -> Result<Self, MatrixError>
pub fn cholesky(&self) -> Result<Self, MatrixError>
Cholesky decomposition: returns lower-triangular L such that A = L·Lᵀ.
§Algorithm
Column-by-column Cholesky-Banachiewicz factorization.
Reference: Golub & Van Loan (1996), Matrix Computations, Algorithm 4.2.1.
§Requirements
Matrix must be symmetric and positive-definite.
§Complexity
O(n³/3).
§Examples
use u_numflow::matrix::Matrix;
let a = Matrix::from_rows(&[
&[4.0, 2.0],
&[2.0, 3.0],
]);
let l = a.cholesky().unwrap();
let llt = l.mul_mat(&l.transpose()).unwrap();
assert!((llt.get(0, 0) - 4.0).abs() < 1e-10);
assert!((llt.get(0, 1) - 2.0).abs() < 1e-10);Sourcepub fn cholesky_solve(&self, b: &[f64]) -> Result<Vec<f64>, MatrixError>
pub fn cholesky_solve(&self, b: &[f64]) -> Result<Vec<f64>, MatrixError>
Solves the linear system A·x = b using Cholesky decomposition.
Equivalent to computing x = A⁻¹·b but more efficient and stable.
§Algorithm
- Decompose A = L·Lᵀ via Cholesky
- Solve L·y = b (forward substitution)
- Solve Lᵀ·x = y (backward substitution)
§Requirements
Matrix must be symmetric positive-definite. b.len() must equal self.rows().
Sourcepub fn eigen_symmetric(&self) -> Result<(Vec<f64>, Matrix), MatrixError>
pub fn eigen_symmetric(&self) -> Result<(Vec<f64>, Matrix), MatrixError>
Eigenvalue decomposition of a real symmetric matrix using the classical Jacobi rotation algorithm.
Returns (eigenvalues, eigenvectors) where eigenvalues are sorted
in descending order and eigenvectors are the corresponding columns
of the returned matrix (column i is the eigenvector for eigenvalue i).
§Algorithm
Cyclic Jacobi rotations zero off-diagonal elements iteratively. Converges quadratically for symmetric matrices.
Reference: Golub & Van Loan (1996), “Matrix Computations”, §8.4
§Complexity
O(n³) per sweep, typically 5–10 sweeps. Best for n < 200.
§Errors
Returns NotSquare if the matrix is not square, or NotSymmetric
if the matrix is not symmetric within tolerance.
§Examples
use u_numflow::matrix::Matrix;
let a = Matrix::from_rows(&[
&[4.0, 1.0],
&[1.0, 3.0],
]);
let (eigenvalues, eigenvectors) = a.eigen_symmetric().unwrap();
// Eigenvalues of [[4,1],[1,3]] are (7+√5)/2 ≈ 4.618 and (7-√5)/2 ≈ 2.382
assert!((eigenvalues[0] - 4.618).abs() < 0.01);
assert!((eigenvalues[1] - 2.382).abs() < 0.01);
// Eigenvectors are orthonormal
let dot: f64 = (0..2).map(|i| eigenvectors.get(i, 0) * eigenvectors.get(i, 1)).sum();
assert!(dot.abs() < 1e-10);