Scilab Website | Contribute with GitLab | Mailing list archives | ATOMS toolboxes
Scilab Online Help
2024.1.0 - Français


taucs_chfact

cholesky factorization of a sparse s.p.d. matrix

Syntax

C_ptr = taucs_chfact(A)

Arguments

A

a sparse real symmetric positive definite (s.p.d.) matrix

C_ptr

a pointer to the Cholesky factors (C,p : A(p,p)=CC')

Description

This function computes a Cholesky factorization of the sparse symmetric positive definite (s.p.d.) matrix A and retrieves at the scilab level, a pointer (C_ptr) to a handle of the Cholesky factors (C,p) (the memory used for them is "outside" scilab space).

If your matrix is s.p.d. this function must be used in place of umf_lufact or in place of the scilab function chfact for a gain in speed.

When such a factorization have been computed, a linear system must be solved with taucs_chsolve. To free the memory used by the Cholesky factors, use taucs_chdel(C_ptr) ;to retrieve the Cholesky factors at the scilab level (for example to display their sparse patterns), use taucs_chget; to get some information (number of non zeros in C), use taucs_chinfo. To compute an approximation of the condition number in norm 2 use cond2sp.

Remarks

  • taucs_chfact works only with the upper triangle of A, and the matrix A must be provided either in its complete form (that is with the lower triangle also) or only with its upper triangle;

  • currently taucs_chfact uses the genmmd (generalized minimum degree) algorithm of Liu to find in a first step the permutation p (so as to minimize the fill-in in the factorization); future versions will let the user choose his/her own reordering by providing a supplementary argument p.

Examples

// Example #1 : a small linear test system
// whom solution must be [1;2;3;4;5]
A = sparse( [ 2 -1  0  0  0;
             -1  2 -1  0  0;
              0 -1  2 -1  0;
              0  0 -1  2 -1;
              0  0  0 -1  2] );
b = [0 ; 0; 0; 0; 6];
Cp = taucs_chfact(A);
x = taucs_chsolve(Cp,b)
// do not forget to clear memory with
taucs_chdel(Cp)
// Example #2 a real example
// first load a sparse matrix
[A] = ReadHBSparse(SCI+"/modules/umfpack/demos/bcsstk24.rsa");
// compute the factorization
Cp = taucs_chfact(A);
b = rand(size(A,1),1); // a random rhs
// use taucs_chsolve for solving Ax=b
x = taucs_chsolve(Cp,b);
norm(A*x - b)
// the same with one iterative refinement step
x = taucs_chsolve(Cp,b,A);
norm(A*x - b)
// do not forget to clear memory
taucs_chdel(Cp)

See also

  • taucs_chsolve — solves a linear s.p.d. system A*X = B from Cholesky factors of the sparse A
  • taucs_chdel — utility function used with taucs_chfact
  • taucs_chinfo — get information on Cholesky factors
  • taucs_chget — retrieve the Cholesky factorization at the scilab level
  • cond2sp — computes an approximation of the 2-norm condition number of a s.p.d. sparse matrix
Report an issue
<< taucs_chdel Interface avec UMFPACK (sparse) taucs_chget >>

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:
Mon Jun 17 17:52:28 CEST 2024