Struct russell_sparse::LinSolver
source · pub struct LinSolver<'a> {
pub actual: Box<dyn Send + LinSolTrait + 'a>,
}Expand description
Unifies the access to linear system solvers
Fields§
§actual: Box<dyn Send + LinSolTrait + 'a>Holds the actual implementation
Implementations§
source§impl<'a> LinSolver<'a>
impl<'a> LinSolver<'a>
sourcepub fn compute(
genie: Genie,
x: &mut Vector,
mat: &mut SparseMatrix,
rhs: &Vector,
params: Option<LinSolParams>
) -> Result<Self, StrError>
pub fn compute( genie: Genie, x: &mut Vector, mat: &mut SparseMatrix, rhs: &Vector, params: Option<LinSolParams> ) -> Result<Self, StrError>
Computes the solution of a linear system
Solves the linear system:
A · x = rhs
(m,n) (n) (m)
§Output
x– the vector of unknown values with dimension equal to mat.ncol
§Input
genie– the actual implementation that does all the magicmat– the matrix representing the sparse coefficient matrix A (see Notes below)rhs– the right-hand side vector with know values an dimension equal to coo.nrowverbose– shows messages
§Notes
- For symmetric matrices,
MUMPSrequires that the symmetry/storage be Lower or Full. - For symmetric matrices,
UMFPACKrequires that the symmetry/storage be Full. - This function calls the actual implementation (genie) via the functions
factorize, andsolve. - This function is best for a single-use, whereas the actual solver should be considered for a recurrent use (e.g., inside a loop).
§Examples
use russell_lab::{vec_approx_eq, Vector};
use russell_sparse::prelude::*;
use russell_sparse::StrError;
fn main() -> Result<(), StrError> {
// constants
let ndim = 3; // number of rows = number of columns
let nnz = 5; // number of non-zero values
// allocate the coefficient matrix
let mut mat = SparseMatrix::new_coo(ndim, ndim, nnz, Sym::No)?;
mat.put(0, 0, 0.2)?;
mat.put(0, 1, 0.2)?;
mat.put(1, 0, 0.5)?;
mat.put(1, 1, -0.25)?;
mat.put(2, 2, 0.25)?;
// print matrix
let mut a = mat.as_dense();
let correct = "┌ ┐\n\
│ 0.2 0.2 0 │\n\
│ 0.5 -0.25 0 │\n\
│ 0 0 0.25 │\n\
└ ┘";
assert_eq!(format!("{}", a), correct);
// allocate the right-hand side vector
let rhs = Vector::from(&[1.0, 1.0, 1.0]);
// calculate the solution
let mut x = Vector::new(ndim);
LinSolver::compute(Genie::Umfpack, &mut x, &mut mat, &rhs, None)?;
let correct = vec![3.0, 2.0, 4.0];
vec_approx_eq(&x, &correct, 1e-14);
Ok(())
}Auto Trait Implementations§
impl<'a> Freeze for LinSolver<'a>
impl<'a> !RefUnwindSafe for LinSolver<'a>
impl<'a> Send for LinSolver<'a>
impl<'a> !Sync for LinSolver<'a>
impl<'a> Unpin for LinSolver<'a>
impl<'a> !UnwindSafe for LinSolver<'a>
Blanket Implementations§
source§impl<T> BorrowMut<T> for Twhere
T: ?Sized,
impl<T> BorrowMut<T> for Twhere
T: ?Sized,
source§fn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
Mutably borrows from an owned value. Read more