use faer::{
Mat,
linalg::solvers::Solve,
sparse::linalg::solvers::{Qr, SymbolicQr},
sparse::{SparseColMat, Triplet},
};
use std::ops::Mul;
use crate::linalg::{LinAlgError, LinAlgResult, SparseLinearSolver};
#[derive(Debug, Clone)]
pub struct SparseQRSolver {
factorizer: Option<Qr<usize, f64>>,
symbolic_factorization: Option<SymbolicQr<usize>>,
hessian: Option<SparseColMat<usize, f64>>,
gradient: Option<Mat<f64>>,
covariance_matrix: Option<Mat<f64>>,
standard_errors: Option<Mat<f64>>,
}
impl SparseQRSolver {
pub fn new() -> Self {
SparseQRSolver {
factorizer: None,
symbolic_factorization: None,
hessian: None,
gradient: None,
covariance_matrix: None,
standard_errors: None,
}
}
pub fn hessian(&self) -> Option<&SparseColMat<usize, f64>> {
self.hessian.as_ref()
}
pub fn gradient(&self) -> Option<&Mat<f64>> {
self.gradient.as_ref()
}
pub fn compute_standard_errors(&mut self) -> Option<&Mat<f64>> {
if self.covariance_matrix.is_none() {
self.compute_covariance_matrix();
}
let hessian = self.hessian.as_ref()?;
let n = hessian.ncols();
if let Some(cov) = &self.covariance_matrix {
let mut std_errors = Mat::zeros(n, 1);
for i in 0..n {
let diag_val = cov[(i, i)];
if diag_val >= 0.0 {
std_errors[(i, 0)] = diag_val.sqrt();
} else {
return None;
}
}
self.standard_errors = Some(std_errors);
}
self.standard_errors.as_ref()
}
pub fn reset_covariance(&mut self) {
self.covariance_matrix = None;
self.standard_errors = None;
}
}
impl Default for SparseQRSolver {
fn default() -> Self {
Self::new()
}
}
impl SparseLinearSolver for SparseQRSolver {
fn solve_normal_equation(
&mut self,
residuals: &Mat<f64>,
jacobians: &SparseColMat<usize, f64>,
) -> LinAlgResult<Mat<f64>> {
let jt = jacobians.as_ref().transpose();
let hessian = jt
.to_col_major()
.map_err(|e| {
LinAlgError::MatrixConversion(
"Failed to convert transposed Jacobian to column-major format".to_string(),
)
.log_with_source(e)
})?
.mul(jacobians.as_ref());
let gradient = jacobians.as_ref().transpose().mul(residuals);
let sym = if let Some(ref cached_sym) = self.symbolic_factorization {
cached_sym.clone()
} else {
let new_sym = SymbolicQr::try_new(hessian.symbolic()).map_err(|e| {
LinAlgError::FactorizationFailed("Symbolic QR decomposition failed".to_string())
.log_with_source(e)
})?;
self.symbolic_factorization = Some(new_sym.clone());
new_sym
};
let qr = Qr::try_new_with_symbolic(sym, hessian.as_ref()).map_err(|e| {
LinAlgError::SingularMatrix(
"QR factorization failed (matrix may be singular)".to_string(),
)
.log_with_source(e)
})?;
let dx = qr.solve(-&gradient);
self.hessian = Some(hessian);
self.gradient = Some(gradient);
self.factorizer = Some(qr);
Ok(dx)
}
fn solve_augmented_equation(
&mut self,
residuals: &Mat<f64>,
jacobians: &SparseColMat<usize, f64>,
lambda: f64,
) -> LinAlgResult<Mat<f64>> {
let n = jacobians.ncols();
let jt = jacobians.as_ref().transpose();
let hessian = jt
.to_col_major()
.map_err(|e| {
LinAlgError::MatrixConversion(
"Failed to convert transposed Jacobian to column-major format".to_string(),
)
.log_with_source(e)
})?
.mul(jacobians.as_ref());
let gradient = jacobians.as_ref().transpose().mul(residuals);
let mut lambda_i_triplets = Vec::with_capacity(n);
for i in 0..n {
lambda_i_triplets.push(Triplet::new(i, i, lambda));
}
let lambda_i =
SparseColMat::try_new_from_triplets(n, n, &lambda_i_triplets).map_err(|e| {
LinAlgError::SparseMatrixCreation("Failed to create lambda*I matrix".to_string())
.log_with_source(e)
})?;
let augmented_hessian = hessian.as_ref() + lambda_i;
let sym = if let Some(ref cached_sym) = self.symbolic_factorization {
cached_sym.clone()
} else {
let new_sym = SymbolicQr::try_new(augmented_hessian.symbolic()).map_err(|e| {
LinAlgError::FactorizationFailed(
"Symbolic QR decomposition failed for augmented system".to_string(),
)
.log_with_source(e)
})?;
self.symbolic_factorization = Some(new_sym.clone());
new_sym
};
let qr = Qr::try_new_with_symbolic(sym, augmented_hessian.as_ref()).map_err(|e| {
LinAlgError::SingularMatrix(
"QR factorization failed (matrix may be singular)".to_string(),
)
.log_with_source(e)
})?;
let dx = qr.solve(-&gradient);
self.hessian = Some(hessian);
self.gradient = Some(gradient);
self.factorizer = Some(qr);
Ok(dx)
}
fn get_hessian(&self) -> Option<&SparseColMat<usize, f64>> {
self.hessian.as_ref()
}
fn get_gradient(&self) -> Option<&Mat<f64>> {
self.gradient.as_ref()
}
fn compute_covariance_matrix(&mut self) -> Option<&Mat<f64>> {
if self.factorizer.is_some()
&& self.hessian.is_some()
&& self.covariance_matrix.is_none()
&& let (Some(factorizer), Some(hessian)) = (&self.factorizer, &self.hessian)
{
let n = hessian.ncols();
let identity = Mat::identity(n, n);
let cov_matrix = factorizer.solve(&identity);
self.covariance_matrix = Some(cov_matrix);
}
self.covariance_matrix.as_ref()
}
fn get_covariance_matrix(&self) -> Option<&Mat<f64>> {
self.covariance_matrix.as_ref()
}
}
#[cfg(test)]
mod tests {
use super::*;
const TOLERANCE: f64 = 1e-10;
type TestResult = Result<(), Box<dyn std::error::Error>>;
fn create_test_data()
-> Result<(SparseColMat<usize, f64>, Mat<f64>), faer::sparse::CreationError> {
let triplets = vec![
Triplet::new(0, 0, 1.0),
Triplet::new(0, 1, 0.0),
Triplet::new(0, 2, 1.0),
Triplet::new(1, 0, 0.0),
Triplet::new(1, 1, 1.0),
Triplet::new(1, 2, 1.0),
Triplet::new(2, 0, 1.0),
Triplet::new(2, 1, 1.0),
Triplet::new(2, 2, 0.0),
Triplet::new(3, 0, 1.0),
Triplet::new(3, 1, 0.0),
Triplet::new(3, 2, 0.0),
];
let jacobian = SparseColMat::try_new_from_triplets(4, 3, &triplets)?;
let residuals = Mat::from_fn(4, 1, |i, _| (i + 1) as f64);
Ok((jacobian, residuals))
}
#[test]
fn test_qr_solver_creation() {
let solver = SparseQRSolver::new();
assert!(solver.factorizer.is_none());
let default_solver = SparseQRSolver::default();
assert!(default_solver.factorizer.is_none());
}
#[test]
fn test_qr_solve_normal_equation() -> TestResult {
let mut solver = SparseQRSolver::new();
let (jacobian, residuals) = create_test_data()?;
let solution = solver.solve_normal_equation(&residuals, &jacobian)?;
assert_eq!(solution.nrows(), 3); assert_eq!(solution.ncols(), 1);
assert!(solver.factorizer.is_some());
Ok(())
}
#[test]
fn test_qr_factorizer_caching() -> TestResult {
let mut solver = SparseQRSolver::new();
let (jacobian, residuals) = create_test_data()?;
let sol1 = solver.solve_normal_equation(&residuals, &jacobian)?;
assert!(solver.factorizer.is_some());
let sol2 = solver.solve_normal_equation(&residuals, &jacobian)?;
for i in 0..sol1.nrows() {
assert!((sol1[(i, 0)] - sol2[(i, 0)]).abs() < TOLERANCE);
}
Ok(())
}
#[test]
fn test_qr_solve_augmented_equation() -> TestResult {
let mut solver = SparseQRSolver::new();
let (jacobian, residuals) = create_test_data()?;
let lambda = 0.1;
let solution = solver.solve_augmented_equation(&residuals, &jacobian, lambda)?;
assert_eq!(solution.nrows(), 3); assert_eq!(solution.ncols(), 1);
Ok(())
}
#[test]
fn test_qr_augmented_different_lambdas() -> TestResult {
let mut solver = SparseQRSolver::new();
let (jacobian, residuals) = create_test_data()?;
let lambda1 = 0.01;
let lambda2 = 1.0;
let sol1 = solver.solve_augmented_equation(&residuals, &jacobian, lambda1)?;
let sol2 = solver.solve_augmented_equation(&residuals, &jacobian, lambda2)?;
let mut different = false;
for i in 0..sol1.nrows() {
if (sol1[(i, 0)] - sol2[(i, 0)]).abs() > TOLERANCE {
different = true;
break;
}
}
assert!(
different,
"Solutions should differ with different lambda values"
);
Ok(())
}
#[test]
fn test_qr_rank_deficient_matrix() -> TestResult {
let mut solver = SparseQRSolver::new();
let triplets = vec![
Triplet::new(0, 0, 1.0),
Triplet::new(0, 1, 2.0),
Triplet::new(0, 2, 3.0),
Triplet::new(1, 0, 2.0),
Triplet::new(1, 1, 4.0),
Triplet::new(1, 2, 6.0), Triplet::new(2, 0, 0.0),
Triplet::new(2, 1, 0.0),
Triplet::new(2, 2, 1.0),
];
let jacobian = SparseColMat::try_new_from_triplets(3, 3, &triplets)?;
let residuals = Mat::from_fn(3, 1, |i, _| i as f64);
let result = solver.solve_normal_equation(&residuals, &jacobian);
assert!(result.is_ok());
Ok(())
}
#[test]
fn test_qr_augmented_system_structure() -> TestResult {
let mut solver = SparseQRSolver::new();
let triplets = vec![
Triplet::new(0, 0, 1.0),
Triplet::new(0, 1, 0.0),
Triplet::new(1, 0, 0.0),
Triplet::new(1, 1, 1.0),
];
let jacobian = SparseColMat::try_new_from_triplets(2, 2, &triplets)?;
let residuals = Mat::from_fn(2, 1, |i, _| (i + 1) as f64);
let lambda = 0.5;
let solution = solver.solve_augmented_equation(&residuals, &jacobian, lambda)?;
assert_eq!(solution.nrows(), 2); assert_eq!(solution.ncols(), 1);
Ok(())
}
#[test]
fn test_qr_numerical_accuracy() -> TestResult {
let mut solver = SparseQRSolver::new();
let triplets = vec![
Triplet::new(0, 0, 1.0),
Triplet::new(1, 1, 1.0),
Triplet::new(2, 2, 1.0),
];
let jacobian = SparseColMat::try_new_from_triplets(3, 3, &triplets)?;
let residuals = Mat::from_fn(3, 1, |i, _| -((i + 1) as f64));
let solution = solver.solve_normal_equation(&residuals, &jacobian)?;
for i in 0..3 {
let expected = (i + 1) as f64;
assert!(
(solution[(i, 0)] - expected).abs() < TOLERANCE,
"Expected {}, got {}",
expected,
solution[(i, 0)]
);
}
Ok(())
}
#[test]
fn test_qr_solver_clone() {
let solver1 = SparseQRSolver::new();
let solver2 = solver1.clone();
assert!(solver1.factorizer.is_none());
assert!(solver2.factorizer.is_none());
}
#[test]
fn test_qr_zero_lambda_augmented() -> TestResult {
let mut solver = SparseQRSolver::new();
let (jacobian, residuals) = create_test_data()?;
let normal_sol = solver.solve_normal_equation(&residuals, &jacobian)?;
let augmented_sol = solver.solve_augmented_equation(&residuals, &jacobian, 0.0)?;
for i in 0..normal_sol.nrows() {
assert!(
(normal_sol[(i, 0)] - augmented_sol[(i, 0)]).abs() < 1e-8,
"Zero lambda augmented should match normal equation"
);
}
Ok(())
}
#[test]
fn test_qr_covariance_computation() -> TestResult {
let mut solver = SparseQRSolver::new();
let (jacobian, residuals) = create_test_data()?;
solver.solve_normal_equation(&residuals, &jacobian)?;
let cov_matrix = solver.compute_covariance_matrix();
assert!(cov_matrix.is_some());
if let Some(cov) = cov_matrix {
assert_eq!(cov.nrows(), 3); assert_eq!(cov.ncols(), 3);
for i in 0..3 {
for j in 0..3 {
assert!(
(cov[(i, j)] - cov[(j, i)]).abs() < TOLERANCE,
"Covariance matrix should be symmetric"
);
}
}
for i in 0..3 {
assert!(
cov[(i, i)] > 0.0,
"Diagonal elements (variances) should be positive"
);
}
}
Ok(())
}
#[test]
fn test_qr_standard_errors_computation() -> TestResult {
let mut solver = SparseQRSolver::new();
let (jacobian, residuals) = create_test_data()?;
solver.solve_normal_equation(&residuals, &jacobian)?;
solver.compute_standard_errors();
assert!(solver.covariance_matrix.is_some());
assert!(solver.standard_errors.is_some());
if let (Some(cov), Some(errors)) = (&solver.covariance_matrix, &solver.standard_errors) {
assert_eq!(errors.nrows(), 3); assert_eq!(errors.ncols(), 1);
for i in 0..3 {
assert!(errors[(i, 0)] > 0.0, "Standard errors should be positive");
}
for i in 0..3 {
let expected_std_error = cov[(i, i)].sqrt();
assert!(
(errors[(i, 0)] - expected_std_error).abs() < TOLERANCE,
"Standard error should equal sqrt of covariance diagonal"
);
}
}
Ok(())
}
#[test]
fn test_qr_covariance_well_conditioned() -> TestResult {
let mut solver = SparseQRSolver::new();
let triplets = vec![
Triplet::new(0, 0, 2.0),
Triplet::new(0, 1, 0.0),
Triplet::new(1, 0, 0.0),
Triplet::new(1, 1, 3.0),
];
let jacobian = SparseColMat::try_new_from_triplets(2, 2, &triplets)?;
let residuals = Mat::from_fn(2, 1, |i, _| (i + 1) as f64);
solver.solve_normal_equation(&residuals, &jacobian)?;
let cov_matrix = solver.compute_covariance_matrix();
assert!(cov_matrix.is_some());
if let Some(cov) = cov_matrix {
assert!((cov[(0, 0)] - 0.25).abs() < TOLERANCE);
assert!((cov[(1, 1)] - 1.0 / 9.0).abs() < TOLERANCE);
assert!(cov[(0, 1)].abs() < TOLERANCE);
assert!(cov[(1, 0)].abs() < TOLERANCE);
}
Ok(())
}
#[test]
fn test_qr_covariance_caching() -> TestResult {
let mut solver = SparseQRSolver::new();
let (jacobian, residuals) = create_test_data()?;
solver.solve_normal_equation(&residuals, &jacobian)?;
solver.compute_covariance_matrix();
assert!(solver.covariance_matrix.is_some());
if let Some(cov1) = &solver.covariance_matrix {
let cov1_ptr = cov1.as_ptr();
solver.compute_covariance_matrix();
assert!(solver.covariance_matrix.is_some());
if let Some(cov2) = &solver.covariance_matrix {
let cov2_ptr = cov2.as_ptr();
assert_eq!(cov1_ptr, cov2_ptr, "Covariance matrix should be cached");
}
}
Ok(())
}
#[test]
fn test_qr_covariance_singular_system() -> TestResult {
let mut solver = SparseQRSolver::new();
let triplets = vec![
Triplet::new(0, 0, 1.0),
Triplet::new(0, 1, 2.0),
Triplet::new(1, 0, 2.0),
Triplet::new(1, 1, 4.0), ];
let jacobian = SparseColMat::try_new_from_triplets(2, 2, &triplets)?;
let residuals = Mat::from_fn(2, 1, |i, _| i as f64);
let result = solver.solve_normal_equation(&residuals, &jacobian);
if result.is_ok() {
let cov_matrix = solver.compute_covariance_matrix();
if let Some(cov) = cov_matrix {
assert!(cov.nrows() == 2);
assert!(cov.ncols() == 2);
}
}
Ok(())
}
}