1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538
|
<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "https://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">
<html xmlns="http://www.w3.org/1999/xhtml">
<head>
<meta http-equiv="Content-Type" content="text/xhtml;charset=UTF-8"/>
<meta http-equiv="X-UA-Compatible" content="IE=9"/>
<meta name="generator" content="Doxygen 1.8.15"/>
<meta name="viewport" content="width=device-width, initial-scale=1"/>
<title>SuperLU: SRC/sgssvx.c File Reference</title>
<link href="tabs.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="jquery.js"></script>
<script type="text/javascript" src="dynsections.js"></script>
<link href="doxygen.css" rel="stylesheet" type="text/css" />
</head>
<body>
<div id="top"><!-- do not remove this div, it is closed by doxygen! -->
<div id="titlearea">
<table cellspacing="0" cellpadding="0">
<tbody>
<tr style="height: 56px;">
<td id="projectalign" style="padding-left: 0.5em;">
<div id="projectname">SuperLU
 <span id="projectnumber">5.3.0</span>
</div>
</td>
</tr>
</tbody>
</table>
</div>
<!-- end header part -->
<!-- Generated by Doxygen 1.8.15 -->
<script type="text/javascript" src="menudata.js"></script>
<script type="text/javascript" src="menu.js"></script>
<script type="text/javascript">
/* @license magnet:?xt=urn:btih:cf05388f2679ee054f2beb29a391d25f4e673ac3&dn=gpl-2.0.txt GPL-v2 */
$(function() {
initMenu('',false,false,'search.php','Search');
});
/* @license-end */</script>
<div id="main-nav"></div>
<div id="nav-path" class="navpath">
<ul>
<li class="navelem"><a class="el" href="dir_1e771ff450ae847412a8c28572c155bb.html">SRC</a></li> </ul>
</div>
</div><!-- top -->
<div class="header">
<div class="summary">
<a href="#func-members">Functions</a> </div>
<div class="headertitle">
<div class="title">sgssvx.c File Reference</div> </div>
</div><!--header-->
<div class="contents">
<p>Solves the system of linear equations A*X=B or A'*X=B.
<a href="#details">More...</a></p>
<div class="textblock"><code>#include "<a class="el" href="slu__sdefs_8h_source.html">slu_sdefs.h</a>"</code><br />
</div><div class="textblock"><div class="dynheader">
Include dependency graph for sgssvx.c:</div>
<div class="dyncontent">
<div class="center"><img src="sgssvx_8c__incl.png" border="0" usemap="#SRC_2sgssvx_8c" alt=""/></div>
<map name="SRC_2sgssvx_8c" id="SRC_2sgssvx_8c">
<area shape="rect" title="Solves the system of linear equations A*X=B or A'*X=B." alt="" coords="267,5,372,32"/>
<area shape="rect" href="slu__sdefs_8h.html" title="Header file for real operations." alt="" coords="275,80,364,107"/>
<area shape="rect" title=" " alt="" coords="5,155,68,181"/>
<area shape="rect" title=" " alt="" coords="93,155,157,181"/>
<area shape="rect" title=" " alt="" coords="95,229,157,256"/>
<area shape="rect" title=" " alt="" coords="181,229,245,256"/>
<area shape="rect" title=" " alt="" coords="432,155,497,181"/>
<area shape="rect" title=" " alt="" coords="552,229,617,256"/>
<area shape="rect" href="slu__Cnames_8h.html" title="Macros defining how C routines will be called." alt="" coords="522,155,629,181"/>
<area shape="rect" href="supermatrix_8h.html" title="Matrix type definitions." alt="" coords="653,155,756,181"/>
<area shape="rect" href="slu__util_8h.html" title="Utility header file." alt="" coords="282,155,357,181"/>
<area shape="rect" title=" " alt="" coords="458,229,527,256"/>
<area shape="rect" href="superlu__enum__consts_8h.html" title="enum constants header file" alt="" coords="269,229,433,256"/>
</map>
</div>
</div><table class="memberdecls">
<tr class="heading"><td colspan="2"><h2 class="groupheader"><a name="func-members"></a>
Functions</h2></td></tr>
<tr class="memitem:af9f64ed49b9ab2d93af6a1100e38ae6d"><td class="memItemLeft" align="right" valign="top">void </td><td class="memItemRight" valign="bottom"><a class="el" href="sgssvx_8c.html#af9f64ed49b9ab2d93af6a1100e38ae6d">sgssvx</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#ac900805a486cbb8489e3c176ed6e0d8e">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, float *ferr, float *berr, <a class="el" href="structGlobalLU__t.html">GlobalLU_t</a> *Glu, <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>
<tr class="separator:af9f64ed49b9ab2d93af6a1100e38ae6d"><td class="memSeparator" colspan="2"> </td></tr>
</table>
<a name="details" id="details"></a><h2 class="groupheader">Detailed Description</h2>
<div class="textblock"><p>Copyright (c) 2003, The Regents of the University of California, through Lawrence Berkeley National Laboratory (subject to receipt of any required approvals from U.S. Dept. of Energy)</p>
<p>All rights reserved.</p>
<p>The source code is distributed under BSD license, see the file License.txt at the top-level directory.</p>
<pre>
-- SuperLU routine (version 3.0) --
Univ. of California Berkeley, Xerox Palo Alto Research Center,
and Lawrence Berkeley National Lab.
October 15, 2003
</pre> </div><h2 class="groupheader">Function Documentation</h2>
<a id="af9f64ed49b9ab2d93af6a1100e38ae6d"></a>
<h2 class="memtitle"><span class="permalink"><a href="#af9f64ed49b9ab2d93af6a1100e38ae6d">◆ </a></span>sgssvx()</h2>
<div class="memitem">
<div class="memproto">
<table class="memname">
<tr>
<td class="memname">void sgssvx </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">float * </td>
<td class="paramname"><em>ferr</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">float * </td>
<td class="paramname"><em>berr</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype"><a class="el" href="structGlobalLU__t.html">GlobalLU_t</a> * </td>
<td class="paramname"><em>Glu</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>
</tr>
<tr>
<td></td>
<td>)</td>
<td></td><td></td>
</tr>
</table>
</div><div class="memdoc">
<pre>
Purpose
=======</pre><pre>SGSSVX solves the system of linear equations A*X=B or A'*X=B, using
the LU factorization from <a class="el" href="sgstrf_8c.html#a7f9874cec10809f11998cc3d9cb88f8b">sgstrf()</a>. Error bounds on the solution and
a condition estimate are also provided. It performs the following steps:</pre><pre> 1. If A is stored column-wise (A->Stype = SLU_NC):</pre><pre> 1.1. If options->Equil = YES, 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><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><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><pre> 1.4. Compute the reciprocal pivot growth factor.</pre><pre> 1.5. If some U(i,i) = 0, so that U is exactly singular, then the
routine returns with info = i. Otherwise, the factored form of
A is used to estimate the condition number of the matrix A. 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 and computes error bounds
as described below.</pre><pre> 1.6. The system of equations is solved for X using the factored form
of A.</pre><pre> 1.7. If options->IterRefine != NOREFINE, iterative refinement is
applied to improve the computed solution matrix and calculate
error bounds and backward error estimates for it.</pre><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><pre> 2. If A is stored row-wise (A->Stype = SLU_NR), apply the above algorithm
to the transpose of A:</pre><pre> 2.1. If options->Equil = YES, 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><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><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><pre> 2.4. Compute the reciprocal pivot growth factor.</pre><pre> 2.5. If some U(i,i) = 0, so that U is exactly singular, then the
routine returns with info = i. Otherwise, the factored form
of transpose(A) is used to estimate the condition number of the
matrix A. If the reciprocal of the condition number
is less than machine precision, info = A->nrow+1 is returned as
a warning, but the routine still goes on to solve for X and
computes error bounds as described below.</pre><pre> 2.6. The system of equations is solved for X using the factored form
of transpose(A).</pre><pre> 2.7. If options->IterRefine != NOREFINE, iterative refinement is
applied to improve the computed solution matrix and calculate
error bounds and backward error estimates for it.</pre><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><pre> See <a class="el" href="supermatrix_8h.html" title="Matrix type definitions.">supermatrix.h</a> for the definition of '<a class="el" href="structSuperMatrix.html">SuperMatrix</a>' structure.</pre><pre>Arguments
=========</pre><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><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_D, Mtype = SLU_GE.
In the future, more general A may be handled.</pre><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.
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><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><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><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><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><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><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><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><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><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><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><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><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><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><pre> See argument 'mem_usage' for memory usage statistics.</pre><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><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><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><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><pre>FERR (output) float*, dimension (B->ncol)
The estimated forward error bound for each solution vector
X(j) (the j-th column of the solution matrix X).
If XTRUE is the true solution corresponding to X(j), FERR(j)
is an estimated upper bound for the magnitude of the largest
element in (X(j) - XTRUE) divided by the magnitude of the
largest element in X(j). The estimate is as reliable as
the estimate for RCOND, and is almost always a slight
overestimate of the true error.
If options->IterRefine = NOREFINE, ferr = 1.0.</pre><pre>BERR (output) float*, dimension (B->ncol)
The componentwise relative backward error of each solution
vector X(j) (i.e., the smallest relative change in
any element of A or B that makes X(j) an exact solution).
If options->IterRefine = NOREFINE, berr = 1.0.</pre><pre>Glu (input/output) <a class="el" href="structGlobalLU__t.html">GlobalLU_t</a> *
If options->Fact == SamePattern_SameRowPerm, it is an input;
The matrix A will be factorized assuming that a
factorization of a matrix with the same sparsity pattern
and similar numerical values was performed prior to this one.
Therefore, this factorization will reuse both row and column
scaling factors R and C, both row and column permutation
vectors perm_r and perm_c, and the L & U data structures
set up from the previous factorization.
Otherwise, it is an output.</pre><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\U 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><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 '<a class="el" href="structSuperLUStat__t.html">SuperLUStat_t</a>'.</pre><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: U(i,i) is exactly zero. The factorization has
been completed, but the factor U is exactly
singular, so the solution and error bounds
could not be computed.
= 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 class="dynheader">
Here is the call graph for this function:</div>
<div class="dyncontent">
<div class="center"><img src="sgssvx_8c_af9f64ed49b9ab2d93af6a1100e38ae6d_cgraph.png" border="0" usemap="#sgssvx_8c_af9f64ed49b9ab2d93af6a1100e38ae6d_cgraph" alt=""/></div>
<map name="sgssvx_8c_af9f64ed49b9ab2d93af6a1100e38ae6d_cgraph" id="sgssvx_8c_af9f64ed49b9ab2d93af6a1100e38ae6d_cgraph">
<area shape="rect" title=" " alt="" coords="5,487,68,513"/>
<area shape="rect" href="slangs_8c.html#a201bfd9f2017cf5904aada9f21f23ab2" title=" " alt="" coords="184,5,244,32"/>
<area shape="rect" href="slu__cdefs_8h.html#a1b49252f1cab66e35ac47ac1afb2adec" title=" " alt="" coords="417,107,479,133"/>
<area shape="rect" href="input__error_8c.html#a5832852b8f6777484b6f55dd79e67d86" title=" " alt="" coords="597,259,684,285"/>
<area shape="rect" href="slu__util_8h.html#a72be96e75e58564c4322ef9ef73ca65f" title=" " alt="" coords="607,689,675,716"/>
<area shape="rect" href="slu__sdefs_8h.html#abb3d30eea43abc536793244e7564e70d" title="Supernodal LU factor related." alt="" coords="360,1677,536,1704"/>
<area shape="rect" href="slu__util_8h.html#a0c6777573bbfe81917cd381e0090d355" title="Timer function." alt="" coords="389,1728,507,1755"/>
<area shape="rect" href="sgsequ_8c.html#ad8a808e807e38c32c08cfbeadb088f08" title="Driver related." alt="" coords="181,208,247,235"/>
<area shape="rect" href="slaqgs_8c.html#af44216962efdebc7e1117b273743e84f" title=" " alt="" coords="184,107,244,133"/>
<area shape="rect" href="get__perm__c_8c.html#aaecb6e6e7a3e97356050bcfdf2573796" title=" " alt="" coords="168,1804,260,1831"/>
<area shape="rect" href="slu__util_8h.html#adf9c573cbfb4520a5ea820702d27cfa5" title=" " alt="" coords="167,1981,261,2008"/>
<area shape="rect" href="sgstrf_8c.html#a7f9874cec10809f11998cc3d9cb88f8b" title=" " alt="" coords="187,1095,241,1121"/>
<area shape="rect" href="slu__sdefs_8h.html#acb8787465a6296109b9a306d5a315ff8" title=" " alt="" coords="163,157,265,184"/>
<area shape="rect" href="sgscon_8c.html#a76b21c7561d5bce81821a76c3465601b" title=" " alt="" coords="416,360,480,387"/>
<area shape="rect" href="sgstrs_8c.html#a9b6e1e555af9cf109ef3a584054a91e2" title=" " alt="" coords="419,309,477,336"/>
<area shape="rect" href="sgsrfs_8c.html#aa619758588187cd5ad69a10a808d18f6" title=" " alt="" coords="185,309,243,336"/>
<area shape="rect" href="slu__sdefs_8h.html#a1357f9a3b2ffb9522883ad84affa63e3" title=" " alt="" coords="163,613,265,640"/>
<area shape="rect" href="slu__util_8h.html#a4de38e1c0ef18dd0791cb206c7f5348f" title="A is of type Stype==NCP." alt="" coords="116,2032,312,2059"/>
<area shape="rect" href="slu__util_8h.html#a2c43be55861c6e4ee5b806ac16cc382c" title="Deallocate the structure pointing to the actual storage of the matrix." alt="" coords="141,2083,287,2125"/>
<area shape="rect" href="SRC_2sp__ienv_8c.html#a89d63af74d9accdbcd7e859b687130cc" title=" " alt="" coords="745,689,832,716"/>
<area shape="rect" href="get__perm__c_8c.html#a90f30e2b284864f6a800a98ceaff8fbc" title=" " alt="" coords="419,1779,477,1805"/>
<area shape="rect" href="get__perm__c_8c.html#a486ee50799ff66abe91efa46a5950a57" title=" " alt="" coords="408,1931,488,1957"/>
<area shape="rect" href="get__perm__c_8c.html#ae92c26cd488b7a86b8277cee2773d8ef" title=" " alt="" coords="403,1829,493,1856"/>
<area shape="rect" href="get__perm__c_8c.html#a792508355b6bef974fcd9e214de40c8e" title=" " alt="" coords="407,1880,489,1907"/>
<area shape="rect" href="slu__util_8h.html#a8a3ba6cbe163f9c12f6f10ee8ba98fc7" title=" " alt="" coords="399,2032,497,2059"/>
<area shape="rect" href="sp__preorder_8c.html#ac79059104ae6abf212c41986820d358c" title=" " alt="" coords="401,2083,495,2109"/>
<area shape="rect" href="sp__coletree_8c.html#a657d6b291654432e815392c2a00d2b84" title=" " alt="" coords="399,2133,497,2160"/>
<area shape="rect" href="slu__util_8h.html#af8198f26bef3c82fbb8601fc5a8e0d9e" title=" " alt="" coords="403,2184,493,2211"/>
<area shape="rect" href="slu__util_8h.html#a44084fde835d2ccaa25e9fd942a72b7a" title=" " alt="" coords="587,1576,694,1603"/>
<area shape="rect" href="slu__sdefs_8h.html#af68715ec86cde90aa31fec07164d6ea6" title="Memory-related." alt="" coords="401,512,495,539"/>
<area shape="rect" href="memory_8c.html#a49bbe20102e5b541c8e8963afa2bd46a" title=" " alt="" coords="603,461,678,488"/>
<area shape="rect" href="memory_8c.html#adbbe5a57b4ed64564c887fb52d798c54" title="Set up pointers for integer working arrays." alt="" coords="409,1171,487,1197"/>
<area shape="rect" href="slu__util_8h.html#ab0dfb6551008bcad5e758defdbd13006" title="Fills an integer array with a given value." alt="" coords="621,1145,660,1172"/>
<area shape="rect" href="slu__sdefs_8h.html#ab5b2859bf1ef1900506dfa702574c6ad" title="Set up pointers for real working arrays." alt="" coords="403,715,493,741"/>
<area shape="rect" href="heap__relax__snode_8c.html#a059d36bb76b7562c9bb2cbd7870e7ffe" title=" " alt="" coords="382,1069,514,1096"/>
<area shape="rect" href="relax__snode_8c.html#ad70bc12cb9031ab8aba4a37a18be46e3" title=" " alt="" coords="401,1221,495,1248"/>
<area shape="rect" href="slu__sdefs_8h.html#ad9d54c8dfc11f1e034b4b7175be60ffb" title=" " alt="" coords="403,917,493,944"/>
<area shape="rect" href="slu__sdefs_8h.html#a9af26d0426eb0bb63755880f2e67e7b7" title="Expand the data structures for L and U during the factorization." alt="" coords="584,917,697,944"/>
<area shape="rect" href="slu__sdefs_8h.html#a60e60255360fae0b1458da070690a3a2" title="Performs numeric block updates within the relaxed snode." alt="" coords="395,1272,501,1299"/>
<area shape="rect" href="slu__sdefs_8h.html#ad7ddf03faedae25b4d73e0b6b33bf50c" title=" " alt="" coords="416,1323,480,1349"/>
<area shape="rect" href="slu__sdefs_8h.html#a297455c494a78c098b2bf418edbc6b16" title="Diagnostic print of column "jcol" in the U/L factor." alt="" coords="399,1373,497,1400"/>
<area shape="rect" href="slu__sdefs_8h.html#a77baf210393e04fa71d4e73b5e60e556" title=" " alt="" coords="405,1424,491,1451"/>
<area shape="rect" href="slu__sdefs_8h.html#a51486936a9ff5079afed80eb5bf8a3e0" title=" " alt="" coords="397,765,499,792"/>
<area shape="rect" href="scolumn__dfs_8c.html#ac9a044320fe8bfbb051a344686a4cb7d" title=" " alt="" coords="399,816,497,843"/>
<area shape="rect" href="scolumn__bmod_8c.html#ac89043410fd16fe2b8d3b2c902fec9f7" title=" " alt="" coords="392,968,504,995"/>
<area shape="rect" href="scopy__to__ucol_8c.html#ab44d465713c602e68295999c003daf7d" title=" " alt="" coords="395,1019,501,1045"/>
<area shape="rect" href="slu__sdefs_8h.html#acf9da2c45289246ef663fc4a96d1ad78" title=" " alt="" coords="413,1475,483,1501"/>
<area shape="rect" href="slu__util_8h.html#ab71db926d60d7b8fd739df197b766366" title="Reset repfnz[] for the current column." alt="" coords="400,1525,496,1552"/>
<area shape="rect" href="slu__cdefs_8h.html#a8086902aa8be3fc7d04c3c82ec3a79dc" title="Count the total number of nonzeros in factors L and U, and in the symmetrically reduced L." alt="" coords="413,1576,483,1603"/>
<area shape="rect" href="slu__cdefs_8h.html#a7061332d759d7e4d73c1b2e5cb0bf2bf" title="Fix up the data storage lsub for L-subscripts. It removes the subscript sets for structural pruning,..." alt="" coords="419,1627,477,1653"/>
</map>
</div>
</div>
</div>
</div><!-- contents -->
<!-- start footer part -->
<hr class="footer"/><address class="footer"><small>
Generated by  <a href="http://www.doxygen.org/index.html">
<img class="footer" src="doxygen.png" alt="doxygen"/>
</a> 1.8.15
</small></address>
</body>
</html>
|