Scilab Website | Contribute with GitLab | Mailing list archives | ATOMS toolboxes
Scilab Online Help
5.4.1 - English

Change language to:
Français - 日本語 - Português - Русский

Please note that the recommended version of Scilab is 2025.0.0. This page might be outdated.
See the recommended documentation of this function

Scilab help >> UMFPACK Interface (sparse) > umf_lufact

umf_lufact

lu factorisation of a sparse matrix

Calling Sequence

LU_ptr = umf_lufact(A)

Arguments

A

a sparse, real or complex, square or rectangular, matrix

LU_ptr

a pointer to umf lu factors (L,U,p,q,R)

Description

This function computes a LU factorisation of the sparse matrix A () and return at the scilab level, a pointer (LU_ptr) to a handle of the LU factors (L,U,p,q,R) (the memory used for them is "outside" scilab stack).

This function must be used in place of umfpack if you have multiple linear systems with the same matrix to solve when the rhs are not known at the same time (for instance A x1 = b1 and A x2 = b2 but b2 depends on x1, etc...).

When such a factorisation have been computed, a linear system must be solved with umf_lusolve (in general x = umf_lusolve(LU_ptr, b) but others options are possible, see umf_lusolve. To free the memory used by the LU factors, use umf_ludel(LU_ptr) (umf_ludel); to retrieve the LU factors at the scilab level (for example to display their sparse patterns), use umf_luget; to get some information (number of non zeros in L and U), use umf_luinfo. To compute an approximation of the condition number use condestsp

Examples

// this is the small linear test system from UMFPACK
// whom solution must be [1;2;3;4;5]
A = sparse( [ 2  3  0  0  0;
              3  0  4  0  6; 
              0 -1 -3  2  0; 
              0  0  1  0  0; 
              0  4  2  0  1] );
b = [8 ; 45; -3; 3; 19];
Lup = umf_lufact(A);
x = umf_lusolve(Lup,b)

// solve now A'x=b
x = umf_lusolve(Lup,b,"A''x=b")
norm(A'*x - b)

// do not forget to clear memory with
umf_ludel(Lup)

// a real (but small)  example
// first load a sparse matrix
[A] = ReadHBSparse(SCI+"/modules/umfpack/examples/arc130.rua");
// compute the factorisation
Lup = umf_lufact(A); 
b = rand(size(A,1),1); // a random rhs
// use umf_lusolve for solving Ax=b
x = umf_lusolve(Lup,b);
norm(A*x - b)

// now the same thing with iterative refiment
x = umf_lusolve(Lup,b,"Ax=b",A);
norm(A*x - b)

// solve now the system A'x=b
x = umf_lusolve(Lup,b,"A''x=b");  // without refinement
norm(A'*x - b)
x = umf_lusolve(Lup,b,"A''x=b",A);  // with refinement
norm(A'*x - b)

// do not forget to clear memory
umf_ludel(Lup)

See Also

  • umfpack — solve sparse linear system
  • umf_luget — retrieve lu factors at the Scilab level
  • umf_lusolve — solve a linear sparse system given the LU factors
  • umf_ludel — utility function used with umf_lufact
  • umf_luinfo — get information on LU factors
  • condestsp — estimate the condition number of a sparse matrix
Report an issue
<< umf_ludel UMFPACK Interface (sparse) umf_luget >>

Copyright (c) 2022-2024 (Dassault Systèmes)
Copyright (c) 2017-2022 (ESI Group)
Copyright (c) 2011-2017 (Scilab Enterprises)
Copyright (c) 1989-2012 (INRIA)
Copyright (c) 1989-2007 (ENPC)
with contributors
Last updated:
Tue Apr 02 17:36:22 CEST 2013