1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249
use crate::matrix::Matrix;
use crate::vector::Vector;
use crate::{to_i32, StrError};
extern "C" {
// Computes the solution to a system of linear equations
// <https://www.netlib.org/lapack/explore-html/d8/d72/dgesv_8f.html>
fn c_dgesv(
n: *const i32,
nrhs: *const i32,
a: *mut f64,
lda: *const i32,
ipiv: *mut i32,
b: *mut f64,
ldb: *const i32,
info: *mut i32,
);
}
/// (dgesv) Solves a general linear system (real numbers)
///
/// For a general matrix `a` (square, symmetric, non-symmetric, dense,
/// sparse), find `x` such that:
///
/// ```text
/// a ⋅ x = b
/// (m,m) (m) (m)
/// ```
///
/// However, the right-hand-side will hold the solution:
///
/// ```text
/// b := a⁻¹⋅b == x
/// ```
///
/// The solution is obtained via LU decomposition using Lapack dgesv routine.
///
/// See also: <https://www.netlib.org/lapack/explore-html/d8/d72/dgesv_8f.html>
///
/// # Note
///
/// 1. The matrix `a` will be modified
/// 2. The right-hand-side `b` will contain the solution `x`
///
/// # Examples
///
/// ```
/// use russell_lab::{solve_lin_sys, Matrix, Vector, StrError};
///
/// fn main() -> Result<(), StrError> {
/// // set matrix and right-hand side
/// let mut a = Matrix::from(&[
/// [1.0, 3.0, -2.0],
/// [3.0, 5.0, 6.0],
/// [2.0, 4.0, 3.0],
/// ]);
/// let mut b = Vector::from(&[5.0, 7.0, 8.0]);
///
/// // solve linear system b := a⁻¹⋅b
/// solve_lin_sys(&mut b, &mut a)?;
///
/// // check
/// let x_correct = "┌ ┐\n\
/// │ -15.000 │\n\
/// │ 8.000 │\n\
/// │ 2.000 │\n\
/// └ ┘";
/// assert_eq!(format!("{:.3}", b), x_correct);
/// Ok(())
/// }
/// ```
pub fn solve_lin_sys(b: &mut Vector, a: &mut Matrix) -> Result<(), StrError> {
let (m, n) = a.dims();
if m != n {
return Err("matrix must be square");
}
if b.dim() != m {
return Err("vector has wrong dimension");
}
if m == 0 {
return Ok(());
}
let mut ipiv = vec![0; m];
let m_i32 = to_i32(m);
let nrhs = 1;
let lda = to_i32(m);
let ldb = lda;
let mut info = 0;
unsafe {
c_dgesv(
&m_i32,
&nrhs,
a.as_mut_data().as_mut_ptr(),
&lda,
ipiv.as_mut_ptr(),
b.as_mut_data().as_mut_ptr(),
&ldb,
&mut info,
)
}
if info < 0 {
println!("LAPACK ERROR (dgesv): Argument #{} had an illegal value", -info);
return Err("LAPACK ERROR (dgesv): An argument had an illegal value");
} else if info > 0 {
println!("LAPACK ERROR (dgesv): U({},{}) is exactly zero", info - 1, info - 1);
return Err("LAPACK ERROR (dgesv): The factorization has been completed, but the factor U is exactly singular");
}
Ok(())
}
////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
#[cfg(test)]
mod tests {
use super::{solve_lin_sys, Matrix, Vector};
use crate::{mat_vec_mul, vec_add, vec_approx_eq, vec_norm, Norm};
#[test]
fn solve_lin_sys_fails_on_non_square() {
let mut a = Matrix::new(2, 3);
let mut b = Vector::new(3);
assert_eq!(solve_lin_sys(&mut b, &mut a), Err("matrix must be square"));
}
#[test]
fn solve_lin_sys_fails_on_wrong_dims() {
let mut a = Matrix::new(2, 2);
let mut b = Vector::new(3);
assert_eq!(solve_lin_sys(&mut b, &mut a), Err("vector has wrong dimension"));
}
#[test]
fn solve_lin_sys_0x0_works() {
let mut a = Matrix::new(0, 0);
let mut b = Vector::new(0);
solve_lin_sys(&mut b, &mut a).unwrap();
assert_eq!(b.dim(), 0);
}
#[test]
fn solve_lin_sys_works() {
#[rustfmt::skip]
let mut a = Matrix::from(&[
[2.0, 1.0, 1.0, 3.0, 2.0],
[1.0, 2.0, 2.0, 1.0, 1.0],
[1.0, 2.0, 9.0, 1.0, 5.0],
[3.0, 1.0, 1.0, 7.0, 1.0],
[2.0, 1.0, 5.0, 1.0, 8.0],
]);
#[rustfmt::skip]
let mut b = Vector::from(&[
-2.0,
4.0,
3.0,
-5.0,
1.0,
]);
solve_lin_sys(&mut b, &mut a).unwrap();
#[rustfmt::skip]
let x_correct = &[
-629.0 / 98.0,
237.0 / 49.0,
-53.0 / 49.0,
62.0 / 49.0,
23.0 / 14.0,
];
vec_approx_eq(&b, x_correct, 1e-13);
}
#[test]
fn solve_lin_sys_1_works() {
// example from https://numericalalgorithmsgroup.github.io/LAPACK_Examples/examples/doc/dgesv_example.html
#[rustfmt::skip]
let mut a = Matrix::from(&[
[ 1.80, 2.88, 2.05, -0.89],
[ 5.25, -2.95, -0.95, -3.80],
[ 1.58, -2.69, -2.90, -1.04],
[-1.11, -0.66, -0.59, 0.80],
]);
#[rustfmt::skip]
let mut b = Vector::from(&[
9.52,
24.35,
0.77,
-6.22,
]);
solve_lin_sys(&mut b, &mut a).unwrap();
#[rustfmt::skip]
let x_correct = &[
1.0,
-1.0,
3.0,
-5.0,
];
vec_approx_eq(&b, x_correct, 1e-14);
}
#[test]
fn solve_lin_sys_2_works() {
const TARGET: f64 = 1234.0;
for m in [0, 1, 5, 7, 12_usize] {
// prepare matrix and rhs
let mut a = Matrix::filled(m, m, 1.0);
let mut b = Vector::filled(m, TARGET);
for i in 0..m {
for j in (i + 1)..m {
a.mul(i, j, -1.0);
}
}
// take copies
let a_copy = a.clone();
let b_copy = b.clone();
// solve linear system: b := a⁻¹⋅b == x
solve_lin_sys(&mut b, &mut a).unwrap();
// compare solution
if m == 0 {
let empty: [f64; 0] = [];
assert_eq!(b.as_data(), &empty);
} else {
let mut x_correct = Vector::new(m);
x_correct[0] = TARGET;
assert_eq!(b.as_data(), x_correct.as_data());
}
// check solution a ⋅ x == rhs (with x == b)
let mut rhs = Vector::new(m);
let mut diff = Vector::new(m);
mat_vec_mul(&mut rhs, 1.0, &a_copy, &b).unwrap();
vec_add(&mut diff, 1.0, &rhs, -1.0, &b_copy).unwrap();
assert_eq!(vec_norm(&diff, Norm::Max), 0.0);
}
}
#[test]
fn solve_lin_sys_singular_handles_error() {
let mut a = Matrix::from(&[
[0.0, 0.0], //
[0.0, 1.0], //
]);
let mut b = Vector::from(&[1.0, 1.0]);
assert_eq!(
solve_lin_sys(&mut b, &mut a).err(),
Some("LAPACK ERROR (dgesv): The factorization has been completed, but the factor U is exactly singular")
);
}
}