<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML 4.01 Transitional//EN">
<html><head><meta http-equiv="Content-Type" content="text/html;charset=UTF-8">
<title>SuperLU: SRC/sgsisx.c File Reference</title>
<link href="doxygen.css" rel="stylesheet" type="text/css">
<link href="tabs.css" rel="stylesheet" type="text/css">
</head><body>
<!-- Generated by Doxygen 1.5.5 -->
<div class="navigation" id="top">
<div class="tabs">
<ul>
<li><a href="index.html"><span>Main Page</span></a></li>
<li><a href="annotated.html"><span>Data Structures</span></a></li>
<li class="current"><a href="files.html"><span>Files</span></a></li>
</ul>
</div>
</div>
<div class="contents">
<h1>SRC/sgsisx.c File Reference</h1>Computes an approximate solutions of linear equations A*X=B or A'*X=B. <a href="#_details">More...</a>
<p>
<code>#include "<a class="el" href="slu__sdefs_8h-source.html">slu_sdefs.h</a>"</code><br>
<table border="0" cellpadding="0" cellspacing="0">
<tr><td></td></tr>
<tr><td colspan="2"><br><h2>Functions</h2></td></tr>
<tr><td class="memItemLeft" nowrap align="right" valign="top">void </td><td class="memItemRight" valign="bottom"><a class="el" href="sgsisx_8c.html#7ef921fcca8189c43499e3e89e7e05ce">sgsisx</a> (<a class="el" href="structsuperlu__options__t.html">superlu_options_t</a> *options, <a class="el" href="structSuperMatrix.html">SuperMatrix</a> *<a class="el" href="ilu__zdrop__row_8c.html#c900805a486cbb8489e3c176ed6e0d8e">A</a>, int *perm_c, int *perm_r, int *etree, char *equed, float *R, float *C, <a class="el" href="structSuperMatrix.html">SuperMatrix</a> *L, <a class="el" href="structSuperMatrix.html">SuperMatrix</a> *U, void *work, int lwork, <a class="el" href="structSuperMatrix.html">SuperMatrix</a> *B, <a class="el" href="structSuperMatrix.html">SuperMatrix</a> *X, float *recip_pivot_growth, float *rcond, <a class="el" href="structmem__usage__t.html">mem_usage_t</a> *mem_usage, <a class="el" href="structSuperLUStat__t.html">SuperLUStat_t</a> *stat, int *info)</td></tr>
</table>
<hr><a name="_details"></a><h2>Detailed Description</h2>
<pre>
-- SuperLU routine (version 4.1) --
Lawrence Berkeley National Laboratory.
November, 2010
</pre> <hr><h2>Function Documentation</h2>
<a class="anchor" name="7ef921fcca8189c43499e3e89e7e05ce"></a><!-- doxytag: member="sgsisx.c::sgsisx" ref="7ef921fcca8189c43499e3e89e7e05ce" args="(superlu_options_t *options, SuperMatrix *A, int *perm_c, int *perm_r, int *etree, char *equed, float *R, float *C, SuperMatrix *L, SuperMatrix *U, void *work, int lwork, SuperMatrix *B, SuperMatrix *X, float *recip_pivot_growth, float *rcond, mem_usage_t *mem_usage, SuperLUStat_t *stat, int *info)" -->
<div class="memitem">
<div class="memproto">
<table class="memname">
<tr>
<td class="memname">void sgsisx </td>
<td>(</td>
<td class="paramtype"><a class="el" href="structsuperlu__options__t.html">superlu_options_t</a> * </td>
<td class="paramname"> <em>options</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype"><a class="el" href="structSuperMatrix.html">SuperMatrix</a> * </td>
<td class="paramname"> <em>A</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">int * </td>
<td class="paramname"> <em>perm_c</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">int * </td>
<td class="paramname"> <em>perm_r</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">int * </td>
<td class="paramname"> <em>etree</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">char * </td>
<td class="paramname"> <em>equed</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">float * </td>
<td class="paramname"> <em>R</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">float * </td>
<td class="paramname"> <em>C</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype"><a class="el" href="structSuperMatrix.html">SuperMatrix</a> * </td>
<td class="paramname"> <em>L</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype"><a class="el" href="structSuperMatrix.html">SuperMatrix</a> * </td>
<td class="paramname"> <em>U</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">void * </td>
<td class="paramname"> <em>work</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">int </td>
<td class="paramname"> <em>lwork</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype"><a class="el" href="structSuperMatrix.html">SuperMatrix</a> * </td>
<td class="paramname"> <em>B</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype"><a class="el" href="structSuperMatrix.html">SuperMatrix</a> * </td>
<td class="paramname"> <em>X</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">float * </td>
<td class="paramname"> <em>recip_pivot_growth</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">float * </td>
<td class="paramname"> <em>rcond</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype"><a class="el" href="structmem__usage__t.html">mem_usage_t</a> * </td>
<td class="paramname"> <em>mem_usage</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype"><a class="el" href="structSuperLUStat__t.html">SuperLUStat_t</a> * </td>
<td class="paramname"> <em>stat</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">int * </td>
<td class="paramname"> <em>info</em></td><td> </td>
</tr>
<tr>
<td></td>
<td>)</td>
<td></td><td></td><td width="100%"></td>
</tr>
</table>
</div>
<div class="memdoc">
<p>
<pre>
Purpose
=======</pre><p>
<pre> SGSISX computes an approximate solutions of linear equations
A*X=B or A'*X=B, using the ILU factorization from <a class="el" href="sgsitrf_8c.html#25788392a605519048cafa995b641fcc">sgsitrf()</a>.
An estimation of the condition number is provided.
The routine performs the following steps:</pre><p>
<pre> 1. If A is stored column-wise (A->Stype = SLU_NC):</pre><p>
<pre> 1.1. If options->Equil = YES or options->RowPerm = LargeDiag, scaling
factors are computed to equilibrate the system:
options->Trans = NOTRANS:
diag(R)*A*diag(C) *inv(diag(C))*X = diag(R)*B
options->Trans = TRANS:
(diag(R)*A*diag(C))**T *inv(diag(R))*X = diag(C)*B
options->Trans = CONJ:
(diag(R)*A*diag(C))**H *inv(diag(R))*X = diag(C)*B
Whether or not the system will be equilibrated depends on the
scaling of the matrix A, but if equilibration is used, A is
overwritten by diag(R)*A*diag(C) and B by diag(R)*B
(if options->Trans=NOTRANS) or diag(C)*B (if options->Trans
= TRANS or CONJ).</pre><p>
<pre> 1.2. Permute columns of A, forming A*Pc, where Pc is a permutation
matrix that usually preserves sparsity.
For more details of this step, see <a class="el" href="sp__preorder_8c.html" title="Permute and performs functions on columns of orginal matrix.">sp_preorder.c</a>.</pre><p>
<pre> 1.3. If options->Fact != FACTORED, the LU decomposition is used to
factor the matrix A (after equilibration if options->Equil = YES)
as Pr*A*Pc = L*U, with Pr determined by partial pivoting.</pre><p>
<pre> 1.4. Compute the reciprocal pivot growth factor.</pre><p>
<pre> 1.5. If some U(i,i) = 0, so that U is exactly singular, then the
routine fills a small number on the diagonal entry, that is
U(i,i) = ||A(:,i)||_oo * options->ILU_FillTol ** (1 - i / n),
and info will be increased by 1. The factored form of A is used
to estimate the condition number of the preconditioner. If the
reciprocal of the condition number is less than machine precision,
info = A->ncol+1 is returned as a warning, but the routine still
goes on to solve for X.</pre><p>
<pre> 1.6. The system of equations is solved for X using the factored form
of A.</pre><p>
<pre> 1.7. options->IterRefine is not used</pre><p>
<pre> 1.8. If equilibration was used, the matrix X is premultiplied by
diag(C) (if options->Trans = NOTRANS) or diag(R)
(if options->Trans = TRANS or CONJ) so that it solves the
original system before equilibration.</pre><p>
<pre> 1.9. options for ILU only
1) If options->RowPerm = LargeDiag, MC64 is used to scale and
permute the matrix to an I-matrix, that is Pr*Dr*A*Dc has
entries of modulus 1 on the diagonal and off-diagonal entries
of modulus at most 1. If MC64 fails, <a class="el" href="dgsequ_8c.html#af22b247cc134fb0ba90285e84ccebb4" title="Driver related.">dgsequ()</a> is used to
equilibrate the system.
( Default: LargeDiag )
2) options->ILU_DropTol = tau is the threshold for dropping.
For L, it is used directly (for the whole row in a supernode);
For U, ||A(:,i)||_oo * tau is used as the threshold
for the i-th column.
If a secondary dropping rule is required, tau will
also be used to compute the second threshold.
( Default: 1e-4 )
3) options->ILU_FillFactor = gamma, used as the initial guess
of memory growth.
If a secondary dropping rule is required, it will also
be used as an upper bound of the memory.
( Default: 10 )
4) options->ILU_DropRule specifies the dropping rule.
Option Meaning
====== ===========
DROP_BASIC: Basic dropping rule, supernodal based ILUTP(tau).
DROP_PROWS: Supernodal based ILUTP(p,tau), p = gamma*nnz(A)/n.
DROP_COLUMN: Variant of ILUTP(p,tau), for j-th column,
p = gamma * nnz(A(:,j)).
DROP_AREA: Variation of ILUTP, for j-th column, use
nnz(F(:,1:j)) / nnz(A(:,1:j)) to control memory.
DROP_DYNAMIC: Modify the threshold tau during factorizaion:
If nnz(L(:,1:j)) / nnz(A(:,1:j)) > gamma
tau_L(j) := MIN(tau_0, tau_L(j-1) * 2);
Otherwise
tau_L(j) := MAX(tau_0, tau_L(j-1) / 2);
tau_U(j) uses the similar rule.
NOTE: the thresholds used by L and U are separate.
DROP_INTERP: Compute the second dropping threshold by
interpolation instead of sorting (default).
In this case, the actual fill ratio is not
guaranteed smaller than gamma.
DROP_PROWS, DROP_COLUMN and DROP_AREA are mutually exclusive.
( Default: DROP_BASIC | DROP_AREA )
5) options->ILU_Norm is the criterion of measuring the magnitude
of a row in a supernode of L. ( Default is INF_NORM )
options->ILU_Norm RowSize(x[1:n])
================= ===============
ONE_NORM ||x||_1 / n
TWO_NORM ||x||_2 / sqrt(n)
INF_NORM max{|x[i]|}
6) options->ILU_MILU specifies the type of MILU's variation.
= SILU: do not perform Modified ILU;
= SMILU_1 (not recommended):
U(i,i) := U(i,i) + sum(dropped entries);
= SMILU_2:
U(i,i) := U(i,i) + SGN(U(i,i)) * sum(dropped entries);
= SMILU_3:
U(i,i) := U(i,i) + SGN(U(i,i)) * sum(|dropped entries|);
NOTE: Even SMILU_1 does not preserve the column sum because of
late dropping.
( Default: SILU )
7) options->ILU_FillTol is used as the perturbation when
encountering zero pivots. If some U(i,i) = 0, so that U is
exactly singular, then
U(i,i) := ||A(:,i)|| * options->ILU_FillTol ** (1 - i / n).
( Default: 1e-2 )</pre><p>
<pre> 2. If A is stored row-wise (A->Stype = SLU_NR), apply the above algorithm
to the transpose of A:</pre><p>
<pre> 2.1. If options->Equil = YES or options->RowPerm = LargeDiag, scaling
factors are computed to equilibrate the system:
options->Trans = NOTRANS:
diag(R)*A*diag(C) *inv(diag(C))*X = diag(R)*B
options->Trans = TRANS:
(diag(R)*A*diag(C))**T *inv(diag(R))*X = diag(C)*B
options->Trans = CONJ:
(diag(R)*A*diag(C))**H *inv(diag(R))*X = diag(C)*B
Whether or not the system will be equilibrated depends on the
scaling of the matrix A, but if equilibration is used, A' is
overwritten by diag(R)*A'*diag(C) and B by diag(R)*B
(if trans='N') or diag(C)*B (if trans = 'T' or 'C').</pre><p>
<pre> 2.2. Permute columns of transpose(A) (rows of A),
forming transpose(A)*Pc, where Pc is a permutation matrix that
usually preserves sparsity.
For more details of this step, see <a class="el" href="sp__preorder_8c.html" title="Permute and performs functions on columns of orginal matrix.">sp_preorder.c</a>.</pre><p>
<pre> 2.3. If options->Fact != FACTORED, the LU decomposition is used to
factor the transpose(A) (after equilibration if
options->Fact = YES) as Pr*transpose(A)*Pc = L*U with the
permutation Pr determined by partial pivoting.</pre><p>
<pre> 2.4. Compute the reciprocal pivot growth factor.</pre><p>
<pre> 2.5. If some U(i,i) = 0, so that U is exactly singular, then the
routine fills a small number on the diagonal entry, that is
U(i,i) = ||A(:,i)||_oo * options->ILU_FillTol ** (1 - i / n).
And info will be increased by 1. The factored form of A is used
to estimate the condition number of the preconditioner. If the
reciprocal of the condition number is less than machine precision,
info = A->ncol+1 is returned as a warning, but the routine still
goes on to solve for X.</pre><p>
<pre> 2.6. The system of equations is solved for X using the factored form
of transpose(A).</pre><p>
<pre> 2.7. If options->IterRefine is not used.</pre><p>
<pre> 2.8. If equilibration was used, the matrix X is premultiplied by
diag(C) (if options->Trans = NOTRANS) or diag(R)
(if options->Trans = TRANS or CONJ) so that it solves the
original system before equilibration.</pre><p>
<pre> See <a class="el" href="supermatrix_8h.html" title="Defines matrix types.">supermatrix.h</a> for the definition of 'SuperMatrix' structure.</pre><p>
<pre> Arguments
=========</pre><p>
<pre> options (input) superlu_options_t*
The structure defines the input parameters to control
how the LU decomposition will be performed and how the
system will be solved.</pre><p>
<pre> A (input/output) SuperMatrix*
Matrix A in A*X=B, of dimension (A->nrow, A->ncol). The number
of the linear equations is A->nrow. Currently, the type of A can be:
Stype = SLU_NC or SLU_NR, Dtype = SLU_S, Mtype = SLU_GE.
In the future, more general A may be handled.</pre><p>
<pre> On entry, If options->Fact = FACTORED and equed is not 'N',
then A must have been equilibrated by the scaling factors in
R and/or C.
On exit, A is not modified
if options->Equil = NO, or
if options->Equil = YES but equed = 'N' on exit, or
if options->RowPerm = NO.</pre><p>
<pre> Otherwise, if options->Equil = YES and equed is not 'N',
A is scaled as follows:
If A->Stype = SLU_NC:
equed = 'R': A := diag(R) * A
equed = 'C': A := A * diag(C)
equed = 'B': A := diag(R) * A * diag(C).
If A->Stype = SLU_NR:
equed = 'R': transpose(A) := diag(R) * transpose(A)
equed = 'C': transpose(A) := transpose(A) * diag(C)
equed = 'B': transpose(A) := diag(R) * transpose(A) * diag(C).</pre><p>
<pre> If options->RowPerm = LargeDiag, MC64 is used to scale and permute
the matrix to an I-matrix, that is A is modified as follows:
P*Dr*A*Dc has entries of modulus 1 on the diagonal and
off-diagonal entries of modulus at most 1. P is a permutation
obtained from MC64.
If MC64 fails, <a class="el" href="sgsequ_8c.html#d8a808e807e38c32c08cfbeadb088f08" title="Driver related.">sgsequ()</a> is used to equilibrate the system,
and A is scaled as above, there is no permutation involved.</pre><p>
<pre> perm_c (input/output) int*
If A->Stype = SLU_NC, Column permutation vector of size A->ncol,
which defines the permutation matrix Pc; perm_c[i] = j means
column i of A is in position j in A*Pc.
On exit, perm_c may be overwritten by the product of the input
perm_c and a permutation that postorders the elimination tree
of Pc'*A'*A*Pc; perm_c is not changed if the elimination tree
is already in postorder.</pre><p>
<pre> If A->Stype = SLU_NR, column permutation vector of size A->nrow,
which describes permutation of columns of transpose(A)
(rows of A) as described above.</pre><p>
<pre> perm_r (input/output) int*
If A->Stype = SLU_NC, row permutation vector of size A->nrow,
which defines the permutation matrix Pr, and is determined
by partial pivoting. perm_r[i] = j means row i of A is in
position j in Pr*A.</pre><p>
<pre> If A->Stype = SLU_NR, permutation vector of size A->ncol, which
determines permutation of rows of transpose(A)
(columns of A) as described above.</pre><p>
<pre> If options->Fact = SamePattern_SameRowPerm, the pivoting routine
will try to use the input perm_r, unless a certain threshold
criterion is violated. In that case, perm_r is overwritten by a
new permutation determined by partial pivoting or diagonal
threshold pivoting.
Otherwise, perm_r is output argument.</pre><p>
<pre> etree (input/output) int*, dimension (A->ncol)
Elimination tree of Pc'*A'*A*Pc.
If options->Fact != FACTORED and options->Fact != DOFACT,
etree is an input argument, otherwise it is an output argument.
Note: etree is a vector of parent pointers for a forest whose
vertices are the integers 0 to A->ncol-1; etree[root]==A->ncol.</pre><p>
<pre> equed (input/output) char*
Specifies the form of equilibration that was done.
= 'N': No equilibration.
= 'R': Row equilibration, i.e., A was premultiplied by diag(R).
= 'C': Column equilibration, i.e., A was postmultiplied by diag(C).
= 'B': Both row and column equilibration, i.e., A was replaced
by diag(R)*A*diag(C).
If options->Fact = FACTORED, equed is an input argument,
otherwise it is an output argument.</pre><p>
<pre> R (input/output) float*, dimension (A->nrow)
The row scale factors for A or transpose(A).
If equed = 'R' or 'B', A (if A->Stype = SLU_NC) or transpose(A)
(if A->Stype = SLU_NR) is multiplied on the left by diag(R).
If equed = 'N' or 'C', R is not accessed.
If options->Fact = FACTORED, R is an input argument,
otherwise, R is output.
If options->zFact = FACTORED and equed = 'R' or 'B', each element
of R must be positive.</pre><p>
<pre> C (input/output) float*, dimension (A->ncol)
The column scale factors for A or transpose(A).
If equed = 'C' or 'B', A (if A->Stype = SLU_NC) or transpose(A)
(if A->Stype = SLU_NR) is multiplied on the right by diag(C).
If equed = 'N' or 'R', C is not accessed.
If options->Fact = FACTORED, C is an input argument,
otherwise, C is output.
If options->Fact = FACTORED and equed = 'C' or 'B', each element
of C must be positive.</pre><p>
<pre> L (output) SuperMatrix*
The factor L from the factorization
Pr*A*Pc=L*U (if A->Stype SLU_= NC) or
Pr*transpose(A)*Pc=L*U (if A->Stype = SLU_NR).
Uses compressed row subscripts storage for supernodes, i.e.,
L has types: Stype = SLU_SC, Dtype = SLU_S, Mtype = SLU_TRLU.</pre><p>
<pre> U (output) SuperMatrix*
The factor U from the factorization
Pr*A*Pc=L*U (if A->Stype = SLU_NC) or
Pr*transpose(A)*Pc=L*U (if A->Stype = SLU_NR).
Uses column-wise storage scheme, i.e., U has types:
Stype = SLU_NC, Dtype = SLU_S, Mtype = SLU_TRU.</pre><p>
<pre> work (workspace/output) void*, size (lwork) (in bytes)
User supplied workspace, should be large enough
to hold data structures for factors L and U.
On exit, if fact is not 'F', L and U point to this array.</pre><p>
<pre> lwork (input) int
Specifies the size of work array in bytes.
= 0: allocate space internally by system malloc;
> 0: use user-supplied work array of length lwork in bytes,
returns error if space runs out.
= -1: the routine guesses the amount of space needed without
performing the factorization, and returns it in
mem_usage->total_needed; no other side effects.</pre><p>
<pre> See argument 'mem_usage' for memory usage statistics.</pre><p>
<pre> B (input/output) SuperMatrix*
B has types: Stype = SLU_DN, Dtype = SLU_S, Mtype = SLU_GE.
On entry, the right hand side matrix.
If B->ncol = 0, only LU decomposition is performed, the triangular
solve is skipped.
On exit,
if equed = 'N', B is not modified; otherwise
if A->Stype = SLU_NC:
if options->Trans = NOTRANS and equed = 'R' or 'B',
B is overwritten by diag(R)*B;
if options->Trans = TRANS or CONJ and equed = 'C' of 'B',
B is overwritten by diag(C)*B;
if A->Stype = SLU_NR:
if options->Trans = NOTRANS and equed = 'C' or 'B',
B is overwritten by diag(C)*B;
if options->Trans = TRANS or CONJ and equed = 'R' of 'B',
B is overwritten by diag(R)*B.</pre><p>
<pre> If options->RowPerm = LargeDiag, MC64 is used to scale and permute
the matrix A to an I-matrix. Then, in addition to the scaling
above, B is further permuted by P*B if options->Trans = NOTRANS,
where P is obtained from MC64.</pre><p>
<pre> X (output) SuperMatrix*
X has types: Stype = SLU_DN, Dtype = SLU_S, Mtype = SLU_GE.
If info = 0 or info = A->ncol+1, X contains the solution matrix
to the original system of equations. Note that A and B are modified
on exit if equed is not 'N', and the solution to the equilibrated
system is inv(diag(C))*X if options->Trans = NOTRANS and
equed = 'C' or 'B', or inv(diag(R))*X if options->Trans = 'T' or 'C'
and equed = 'R' or 'B'.</pre><p>
<pre> recip_pivot_growth (output) float*
The reciprocal pivot growth factor max_j( norm(A_j)/norm(U_j) ).
The infinity norm is used. If recip_pivot_growth is much less
than 1, the stability of the LU factorization could be poor.</pre><p>
<pre> rcond (output) float*
The estimate of the reciprocal condition number of the matrix A
after equilibration (if done). If rcond is less than the machine
precision (in particular, if rcond = 0), the matrix is singular
to working precision. This condition is indicated by a return
code of info > 0.</pre><p>
<pre> mem_usage (output) mem_usage_t*
Record the memory usage statistics, consisting of following fields:<ul>
<li>for_lu (float)
The amount of space used in bytes for L data structures.</li><li>total_needed (float)
The amount of space needed in bytes to perform factorization.</li><li>expansions (int)
The number of memory expansions during the LU factorization.</li></ul>
</pre><p>
<pre> stat (output) SuperLUStat_t*
Record the statistics on runtime and floating-point operation count.
See <a class="el" href="slu__util_8h.html" title="Utility header file.">slu_util.h</a> for the definition of 'SuperLUStat_t'.</pre><p>
<pre> info (output) int*
= 0: successful exit
< 0: if info = -i, the i-th argument had an illegal value
> 0: if info = i, and i is
<= A->ncol: number of zero pivots. They are replaced by small
entries due to options->ILU_FillTol.
= A->ncol+1: U is nonsingular, but RCOND is less than machine
precision, meaning that the matrix is singular to
working precision. Nevertheless, the solution and
error bounds are computed because there are a number
of situations where the computed solution can be more
accurate than the value of RCOND would suggest.
> A->ncol+1: number of bytes allocated when memory allocation
failure occurred, plus A->ncol.
</pre>
</div>
</div><p>
</div>
<hr size="1"><address style="text-align: right;"><small>Generated on Mon Nov 22 10:23:48 2010 for SuperLU by
<a href="http://www.doxygen.org/index.html">
<img src="doxygen.png" alt="doxygen" align="middle" border="0"></a> 1.5.5 </small></address>
</body>
</html>