Scilab Home page | Wiki | Bug tracker | Forge | Mailing list archives | ATOMS | File exchange
Please login or create an account
Change language to: Français - Português - Русский - 日本語
Scilab Help >> UMFPACK Interface (sparse) > umf_luget

umf_luget

retrieve lu factors at the Scilab level

Syntax

[L,U,p,q,Rd] = umf_luget(LU_ptr)

Arguments

LU_ptr

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

L,U

scilab sparse matrix

p,q

column vectors storing the permutations

Rd

vector storing the (row) scaling factors

Description

This function may be used if you want to plot the sparse pattern of the lu factors (or if you code something which use the lu factors). The factorization provided by umfpack is of the form:

P R^(-1) A Q = LU

where P and Q are permutation matrices, R is a diagonal matrix (row scaling), L a lower triangular matrix with a diagonal of 1, and U an upper triangular matrix. The function provides the matrices L and U as Sparse scilab matrices but P and Q are given as permutation vectors p and q (in fact p is the permutation associated to P') and Rd is the vector corresponding to the diagonal of R.

Examples

// this is the test matrix from UMFPACK
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] );
Lup = umf_lufact(A);
[L,U,p,q,R] = umf_luget(Lup);
B = A;
for i=1:5, B(i,:) = B(i,:)/R(i); end // apply the row scaling
B(p,q) - L*U  // must be a (quasi) nul matrix

umf_ludel(Lup) // clear memory
// the same with a complex matrix
A = sparse( [ 2+%i  3+2*%i  0      0    0;
              3-%i  0       4+%i   0    6-3*%i;
              0    -1+%i   -3+6*%i 2-%i 0;
              0     0       1-5*%i 0    0;
              0     4       2-%i   0    1] );
Lup = umf_lufact(A);
[L,U,p,q,R] = umf_luget(Lup);
B = A;
for i=1:5
  B(i,:) = B(i,:)/R(i)
end // apply the row scaling
B(p,q) - L*U  // must be a (quasi) nul matrix

umf_ludel(Lup) // clear memory

See Also

  • umfpack — solve sparse linear system
  • umf_lufact — lu factorization of a sparse matrix
  • 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
Scilab Enterprises
Copyright (c) 2011-2015 (Scilab Enterprises)
Copyright (c) 1989-2012 (INRIA)
Copyright (c) 1989-2007 (ENPC)
with contributors
Last updated:
Wed Jun 15 08:27:38 CEST 2016