File: dgsisx_8c.html

package info (click to toggle)
superlu 5.3.0%2Bdfsg1-2
  • links: PTS, VCS
  • area: main
  • in suites: bookworm
  • size: 11,976 kB
  • sloc: ansic: 59,420; makefile: 405; fortran: 185; csh: 141
file content (573 lines) | stat: -rw-r--r-- 37,329 bytes parent folder | download
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
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
<!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/dgsisx.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
   &#160;<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&amp;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">dgsisx.c File Reference</div>  </div>
</div><!--header-->
<div class="contents">

<p>Computes an approximate solutions of linear equations A*X=B or A'*X=B.  
<a href="#details">More...</a></p>
<div class="textblock"><code>#include &quot;<a class="el" href="slu__ddefs_8h_source.html">slu_ddefs.h</a>&quot;</code><br />
</div><div class="textblock"><div class="dynheader">
Include dependency graph for dgsisx.c:</div>
<div class="dyncontent">
<div class="center"><img src="dgsisx_8c__incl.png" border="0" usemap="#SRC_2dgsisx_8c" alt=""/></div>
<map name="SRC_2dgsisx_8c" id="SRC_2dgsisx_8c">
<area shape="rect"  title="Computes an approximate solutions of linear equations A*X=B or A&#39;*X=B." alt="" coords="268,5,371,32"/>
<area shape="rect"  href="slu__ddefs_8h.html" title="Header file for real operations." alt="" coords="274,80,365,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:a895611a0135f2212f505986a6384f1b9"><td class="memItemLeft" align="right" valign="top">void&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="dgsisx_8c.html#a895611a0135f2212f505986a6384f1b9">dgsisx</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, double *R, double *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, double *recip_pivot_growth, double *rcond, <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:a895611a0135f2212f505986a6384f1b9"><td class="memSeparator" colspan="2">&#160;</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 4.2) --
Lawrence Berkeley National Laboratory.
November, 2010
August, 2011
</pre> </div><h2 class="groupheader">Function Documentation</h2>
<a id="a895611a0135f2212f505986a6384f1b9"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a895611a0135f2212f505986a6384f1b9">&#9670;&nbsp;</a></span>dgsisx()</h2>

<div class="memitem">
<div class="memproto">
      <table class="memname">
        <tr>
          <td class="memname">void dgsisx </td>
          <td>(</td>
          <td class="paramtype"><a class="el" href="structsuperlu__options__t.html">superlu_options_t</a> *&#160;</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> *&#160;</td>
          <td class="paramname"><em>A</em>, </td>
        </tr>
        <tr>
          <td class="paramkey"></td>
          <td></td>
          <td class="paramtype">int *&#160;</td>
          <td class="paramname"><em>perm_c</em>, </td>
        </tr>
        <tr>
          <td class="paramkey"></td>
          <td></td>
          <td class="paramtype">int *&#160;</td>
          <td class="paramname"><em>perm_r</em>, </td>
        </tr>
        <tr>
          <td class="paramkey"></td>
          <td></td>
          <td class="paramtype">int *&#160;</td>
          <td class="paramname"><em>etree</em>, </td>
        </tr>
        <tr>
          <td class="paramkey"></td>
          <td></td>
          <td class="paramtype">char *&#160;</td>
          <td class="paramname"><em>equed</em>, </td>
        </tr>
        <tr>
          <td class="paramkey"></td>
          <td></td>
          <td class="paramtype">double *&#160;</td>
          <td class="paramname"><em>R</em>, </td>
        </tr>
        <tr>
          <td class="paramkey"></td>
          <td></td>
          <td class="paramtype">double *&#160;</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> *&#160;</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> *&#160;</td>
          <td class="paramname"><em>U</em>, </td>
        </tr>
        <tr>
          <td class="paramkey"></td>
          <td></td>
          <td class="paramtype">void *&#160;</td>
          <td class="paramname"><em>work</em>, </td>
        </tr>
        <tr>
          <td class="paramkey"></td>
          <td></td>
          <td class="paramtype">int&#160;</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> *&#160;</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> *&#160;</td>
          <td class="paramname"><em>X</em>, </td>
        </tr>
        <tr>
          <td class="paramkey"></td>
          <td></td>
          <td class="paramtype">double *&#160;</td>
          <td class="paramname"><em>recip_pivot_growth</em>, </td>
        </tr>
        <tr>
          <td class="paramkey"></td>
          <td></td>
          <td class="paramtype">double *&#160;</td>
          <td class="paramname"><em>rcond</em>, </td>
        </tr>
        <tr>
          <td class="paramkey"></td>
          <td></td>
          <td class="paramtype"><a class="el" href="structGlobalLU__t.html">GlobalLU_t</a> *&#160;</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> *&#160;</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> *&#160;</td>
          <td class="paramname"><em>stat</em>, </td>
        </tr>
        <tr>
          <td class="paramkey"></td>
          <td></td>
          <td class="paramtype">int *&#160;</td>
          <td class="paramname"><em>info</em>&#160;</td>
        </tr>
        <tr>
          <td></td>
          <td>)</td>
          <td></td><td></td>
        </tr>
      </table>
</div><div class="memdoc">
<pre>
Purpose
=======</pre><pre>DGSISX computes an approximate solutions of linear equations
A*X=B or A'*X=B, using the ILU factorization from <a class="el" href="dgsitrf_8c.html#ac36e874ce15790ca0b059fab7eed0466">dgsitrf()</a>.
An estimation of the condition number is provided. 
The routine performs the following steps:</pre><pre>  1. If A is stored column-wise (A-&gt;Stype = SLU_NC):</pre><pre>     1.1. If options-&gt;Equil = YES or options-&gt;RowPerm = LargeDiag_MC64, scaling
          factors are computed to equilibrate the system:
          options-&gt;Trans = NOTRANS:
         diag(R)*A*diag(C) *inv(diag(C))*X = diag(R)*B
          options-&gt;Trans = TRANS:
         (diag(R)*A*diag(C))**T *inv(diag(R))*X = diag(C)*B
          options-&gt;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-&gt;Trans=NOTRANS) or diag(C)*B (if options-&gt;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-&gt;Fact != FACTORED, the LU decomposition is used to
          factor the matrix A (after equilibration if options-&gt;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 fills a small number on the diagonal entry, that is
        U(i,i) = ||A(:,i)||_oo * options-&gt;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-&gt;ncol+1 is returned as a warning, but the routine still
          goes on to solve for X.</pre><pre>     1.6. The system of equations is solved for X using the factored form
          of A.</pre><pre>     1.7. options-&gt;IterRefine is not used</pre><pre>     1.8. If equilibration was used, the matrix X is premultiplied by
          diag(C) (if options-&gt;Trans = NOTRANS) or diag(R)
          (if options-&gt;Trans = TRANS or CONJ) so that it solves the
          original system before equilibration.</pre><pre>     1.9. options for ILU only
          1) If options-&gt;RowPerm = LargeDiag_MC64, 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#aaf22b247cc134fb0ba90285e84ccebb4" title="Driver related.">dgsequ()</a> is used to
        equilibrate the system.
             ( Default: LargeDiag_MC64 )
          2) options-&gt;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-&gt;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-&gt;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)) &gt; 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-&gt;ILU_Norm is the criterion of measuring the magnitude
        of a row in a supernode of L. ( Default is INF_NORM )
        options-&gt;ILU_Norm       RowSize(x[1:n])
        =================       ===============
        ONE_NORM                ||x||_1 / n
        TWO_NORM                ||x||_2 / sqrt(n)
        INF_NORM                max{|x[i]|}
          6) options-&gt;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) + <a class="el" href="ilu__zpivotL_8c.html#a95ed41486ca0ed53262e4b8934d4afac">SGN(U(i,i))</a> * sum(dropped entries);
        = SMILU_3:
            U(i,i) := U(i,i) + <a class="el" href="ilu__zpivotL_8c.html#a95ed41486ca0ed53262e4b8934d4afac">SGN(U(i,i))</a> * sum(|dropped entries|);
        NOTE: Even SMILU_1 does not preserve the column sum because of
        late dropping.
             ( Default: SILU )
          7) options-&gt;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-&gt;ILU_FillTol ** (1 - i / n).
             ( Default: 1e-2 )</pre><pre>  2. If A is stored row-wise (A-&gt;Stype = SLU_NR), apply the above algorithm
     to the transpose of A:</pre><pre>     2.1. If options-&gt;Equil = YES or options-&gt;RowPerm = LargeDiag_MC64, scaling
          factors are computed to equilibrate the system:
          options-&gt;Trans = NOTRANS:
         diag(R)*A*diag(C) *inv(diag(C))*X = diag(R)*B
          options-&gt;Trans = TRANS:
         (diag(R)*A*diag(C))**T *inv(diag(R))*X = diag(C)*B
          options-&gt;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-&gt;Fact != FACTORED, the LU decomposition is used to
          factor the transpose(A) (after equilibration if
          options-&gt;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 fills a small number on the diagonal entry, that is
         U(i,i) = ||A(:,i)||_oo * options-&gt;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-&gt;ncol+1 is returned as a warning, but the routine still
          goes on to solve for X.</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-&gt;IterRefine is not used.</pre><pre>     2.8. If equilibration was used, the matrix X is premultiplied by
          diag(C) (if options-&gt;Trans = NOTRANS) or diag(R)
          (if options-&gt;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-&gt;nrow, A-&gt;ncol). The number
        of the linear equations is A-&gt;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-&gt;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-&gt;Equil = NO, or
        if options-&gt;Equil = YES but equed = 'N' on exit, or
        if options-&gt;RowPerm = NO.</pre><pre>        Otherwise, if options-&gt;Equil = YES and equed is not 'N',
        A is scaled as follows:
        If A-&gt;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-&gt;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>        If options-&gt;RowPerm = LargeDiag_MC64, 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="dgsequ_8c.html#aaf22b247cc134fb0ba90285e84ccebb4" title="Driver related.">dgsequ()</a> is used to equilibrate the system,
           and A is scaled as above, but no permutation is involved.
           On exit, A is restored to the orginal row numbering, so
           Dr*A*Dc is returned.</pre><pre>perm_c  (input/output) int*
        If A-&gt;Stype = SLU_NC, Column permutation vector of size A-&gt;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-&gt;Stype = SLU_NR, column permutation vector of size A-&gt;nrow,
        which describes permutation of columns of transpose(A) 
        (rows of A) as described above.</pre><pre>perm_r  (input/output) int*
        If A-&gt;Stype = SLU_NC, row permutation vector of size A-&gt;nrow, 
        which defines the permutation matrix Pr, and is determined
        by MC64 first then followed by partial pivoting.
        perm_r[i] = j means row i of A is in position j in Pr*A.</pre><pre>        If A-&gt;Stype = SLU_NR, permutation vector of size A-&gt;ncol, which
        determines permutation of rows of transpose(A)
        (columns of A) as described above.</pre><pre>        If options-&gt;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-&gt;ncol)
        Elimination tree of Pc'*A'*A*Pc.
        If options-&gt;Fact != FACTORED and options-&gt;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-&gt;ncol-1; etree[root]==A-&gt;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-&gt;Fact = FACTORED, equed is an input argument,
        otherwise it is an output argument.</pre><pre>R          (input/output) double*, dimension (A-&gt;nrow)
        The row scale factors for A or transpose(A).
        If equed = 'R' or 'B', A (if A-&gt;Stype = SLU_NC) or transpose(A)
            (if A-&gt;Stype = SLU_NR) is multiplied on the left by diag(R).
        If equed = 'N' or 'C', R is not accessed.
        If options-&gt;Fact = FACTORED, R is an input argument,
            otherwise, R is output.
        If options-&gt;Fact = FACTORED and equed = 'R' or 'B', each element
            of R must be positive.</pre><pre>C          (input/output) double*, dimension (A-&gt;ncol)
        The column scale factors for A or transpose(A).
        If equed = 'C' or 'B', A (if A-&gt;Stype = SLU_NC) or transpose(A)
            (if A-&gt;Stype = SLU_NR) is multiplied on the right by diag(C).
        If equed = 'N' or 'R', C is not accessed.
        If options-&gt;Fact = FACTORED, C is an input argument,
            otherwise, C is output.
        If options-&gt;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-&gt;Stype SLU_= NC) or
            Pr*transpose(A)*Pc=L*U      (if A-&gt;Stype = SLU_NR).
        Uses compressed row subscripts storage for supernodes, i.e.,
        L has types: Stype = SLU_SC, Dtype = SLU_D, Mtype = SLU_TRLU.</pre><pre>U          (output) SuperMatrix*
        The factor U from the factorization
            Pr*A*Pc=L*U         (if A-&gt;Stype = SLU_NC) or
            Pr*transpose(A)*Pc=L*U      (if A-&gt;Stype = SLU_NR).
        Uses column-wise storage scheme, i.e., U has types:
        Stype = SLU_NC, Dtype = SLU_D, 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;
        &gt; 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-&gt;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_D, Mtype = SLU_GE.
        On entry, the right hand side matrix.
        If B-&gt;ncol = 0, only LU decomposition is performed, the triangular
                   solve is skipped.
        On exit,
           if equed = 'N', B is not modified; otherwise
           if A-&gt;Stype = SLU_NC:
         if options-&gt;Trans = NOTRANS and equed = 'R' or 'B',
            B is overwritten by diag(R)*B;
         if options-&gt;Trans = TRANS or CONJ and equed = 'C' of 'B',
            B is overwritten by diag(C)*B;
           if A-&gt;Stype = SLU_NR:
         if options-&gt;Trans = NOTRANS and equed = 'C' or 'B',
            B is overwritten by diag(C)*B;
         if options-&gt;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_D, Mtype = SLU_GE.
        If info = 0 or info = A-&gt;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-&gt;Trans = NOTRANS and
        equed = 'C' or 'B', or inv(diag(R))*X if options-&gt;Trans = 'T' or 'C'
        and equed = 'R' or 'B'.</pre><pre>recip_pivot_growth (output) double*
        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) double*
        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 &gt; 0.</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
        &lt; 0: if info = -i, the i-th argument had an illegal value
        &gt; 0: if info = i, and i is
        &lt;= A-&gt;ncol: number of zero pivots. They are replaced by small
              entries due to options-&gt;ILU_FillTol.
        = A-&gt;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.
        &gt; A-&gt;ncol+1: number of bytes allocated when memory allocation
              failure occurred, plus A-&gt;ncol.
</pre> <div class="dynheader">
Here is the call graph for this function:</div>
<div class="dyncontent">
<div class="center"><img src="dgsisx_8c_a895611a0135f2212f505986a6384f1b9_cgraph.png" border="0" usemap="#dgsisx_8c_a895611a0135f2212f505986a6384f1b9_cgraph" alt=""/></div>
<map name="dgsisx_8c_a895611a0135f2212f505986a6384f1b9_cgraph" id="dgsisx_8c_a895611a0135f2212f505986a6384f1b9_cgraph">
<area shape="rect"  title=" " alt="" coords="5,563,65,589"/>
<area shape="rect"  href="dlangs_8c.html#a75a53f4464b95c63adad9e1f63f44d1c" title=" " alt="" coords="181,5,242,32"/>
<area shape="rect"  href="dmach_8c.html#a9619845d75710daa0aa4549c2862bae5" title=" " alt="" coords="607,157,670,184"/>
<area shape="rect"  href="input__error_8c.html#a5832852b8f6777484b6f55dd79e67d86" title=" " alt="" coords="402,259,489,285"/>
<area shape="rect"  href="slu__util_8h.html#a72be96e75e58564c4322ef9ef73ca65f" title=" " alt="" coords="605,867,673,893"/>
<area shape="rect"  href="dutil_8c.html#a4a177c54dafbe3640c26caa49eeee1de" title="Supernodal LU factor related." alt="" coords="357,1677,533,1704"/>
<area shape="rect"  href="slu__util_8h.html#a0c6777573bbfe81917cd381e0090d355" title="Timer function." alt="" coords="386,1728,505,1755"/>
<area shape="rect"  href="memory_8c.html#a49bbe20102e5b541c8e8963afa2bd46a" title=" " alt="" coords="601,588,676,615"/>
<area shape="rect"  href="dldperm_8c.html#a5c6a8de5e809f2094735965027adf532" title=" " alt="" coords="176,461,247,488"/>
<area shape="rect"  href="dgsequ_8c.html#aaf22b247cc134fb0ba90285e84ccebb4" title="Driver related." alt="" coords="179,247,244,273"/>
<area shape="rect"  href="dlaqgs_8c.html#a07e1fa4926680eb02069087f0aa26fa1" title=" " alt="" coords="181,119,242,145"/>
<area shape="rect"  href="get__perm__c_8c.html#aaecb6e6e7a3e97356050bcfdf2573796" title=" " alt="" coords="165,1804,257,1831"/>
<area shape="rect"  href="slu__util_8h.html#adf9c573cbfb4520a5ea820702d27cfa5" title=" " alt="" coords="165,1981,258,2008"/>
<area shape="rect"  href="dgsitrf_8c.html#ac36e874ce15790ca0b059fab7eed0466" title=" " alt="" coords="182,1171,241,1197"/>
<area shape="rect"  href="dpivotgrowth_8c.html#a770618182a3841e8d10a26a4eb97418a" title=" " alt="" coords="160,183,263,209"/>
<area shape="rect"  href="dgscon_8c.html#a2c7a4267d306243d3ceb15531522033e" title=" " alt="" coords="179,399,244,425"/>
<area shape="rect"  href="dgstrs_8c.html#a6e3eace519372b7dfcd053e0d3614fc1" title=" " alt="" coords="182,297,241,324"/>
<area shape="rect"  href="dmemory_8c.html#a61aaccf587a78d15d79c4cc79f80e8b0" title=" " alt="" coords="149,715,273,741"/>
<area shape="rect"  href="slu__util_8h.html#a4de38e1c0ef18dd0791cb206c7f5348f" title="A is of type Stype==NCP." alt="" coords="113,2032,309,2059"/>
<area shape="rect"  href="slu__util_8h.html#a2c43be55861c6e4ee5b806ac16cc382c" title="Deallocate the structure pointing to the actual storage of the matrix." alt="" coords="138,2083,285,2125"/>
<area shape="rect"  href="SRC_2sp__ienv_8c.html#a89d63af74d9accdbcd7e859b687130cc" title=" " alt="" coords="744,867,831,893"/>
<area shape="rect"  href="slu__util_8h.html#ade363dcb4babb66fa0e5f51bd2e6e42c" title=" " alt="" coords="393,512,497,539"/>
<area shape="rect"  href="dldperm_8c.html#a1e6fb0c8dd36aef071ef165136ece781" title=" " alt="" coords="409,411,481,437"/>
<area shape="rect"  href="dldperm_8c.html#a4e7ad390efe21b3fb50607aed500777d" title=" " alt="" coords="407,461,484,488"/>
<area shape="rect"  href="get__perm__c_8c.html#a90f30e2b284864f6a800a98ceaff8fbc" title=" " alt="" coords="416,1779,475,1805"/>
<area shape="rect"  href="get__perm__c_8c.html#a486ee50799ff66abe91efa46a5950a57" title=" " alt="" coords="405,1931,485,1957"/>
<area shape="rect"  href="get__perm__c_8c.html#ae92c26cd488b7a86b8277cee2773d8ef" title=" " alt="" coords="400,1829,491,1856"/>
<area shape="rect"  href="get__perm__c_8c.html#a792508355b6bef974fcd9e214de40c8e" title=" " alt="" coords="405,1880,486,1907"/>
<area shape="rect"  href="slu__util_8h.html#a8a3ba6cbe163f9c12f6f10ee8ba98fc7" title=" " alt="" coords="397,1981,494,2008"/>
<area shape="rect"  href="sp__preorder_8c.html#ac79059104ae6abf212c41986820d358c" title=" " alt="" coords="398,2032,493,2059"/>
<area shape="rect"  href="sp__coletree_8c.html#a657d6b291654432e815392c2a00d2b84" title=" " alt="" coords="396,2083,495,2109"/>
<area shape="rect"  href="slu__util_8h.html#af8198f26bef3c82fbb8601fc5a8e0d9e" title=" " alt="" coords="400,2133,491,2160"/>
<area shape="rect"  href="slu__util_8h.html#a44084fde835d2ccaa25e9fd942a72b7a" title=" " alt="" coords="585,1753,692,1780"/>
<area shape="rect"  href="dmemory_8c.html#ae2ca2ac5e9a763fd3f07487343e4522e" title="Allocate storage for the data structures common to all factor routines." alt="" coords="399,765,492,792"/>
<area shape="rect"  href="memory_8c.html#adbbe5a57b4ed64564c887fb52d798c54" title="Set up pointers for integer working arrays." alt="" coords="407,1627,484,1653"/>
<area shape="rect"  href="slu__util_8h.html#ab0dfb6551008bcad5e758defdbd13006" title="Fills an integer array with a given value." alt="" coords="619,1297,658,1324"/>
<area shape="rect"  href="dmemory_8c.html#aaa5359da217b433b43bf6c8e2d29aa45" title="Set up pointers for real working arrays." alt="" coords="400,968,491,995"/>
<area shape="rect"  href="ilu__heap__relax__snode_8c.html#aac1a978dda622cdb58c3c2eaee4b4030" title=" " alt="" coords="369,1323,521,1349"/>
<area shape="rect"  href="ilu__relax__snode_8c.html#ae0e2bbb8507d800766030635a3bd5a7e" title=" " alt="" coords="387,1221,503,1248"/>
<area shape="rect"  href="mark__relax_8c.html#a5e85b0273eec011f0027d8506a20350e" title=" " alt="" coords="401,1373,489,1400"/>
<area shape="rect"  href="ilu__ddrop__row_8c.html#a380317801e05b11930fd1e094db34179" title=" " alt="" coords="393,1424,498,1451"/>
<area shape="rect"  href="ilu__dsnode__dfs_8c.html#a66dbc4626e59d14b4d3458c4eb841829" title=" " alt="" coords="389,1069,501,1096"/>
<area shape="rect"  href="dmemory_8c.html#a9aff5dfe301496ef7c9234789975c043" title="Expand the data structures for L and U during the factorization." alt="" coords="581,1069,696,1096"/>
<area shape="rect"  href="dsnode__bmod_8c.html#a1466b84198911ff34e828a811e70831e" title="Performs numeric block updates within the relaxed snode." alt="" coords="392,1475,499,1501"/>
<area shape="rect"  href="ilu__dpivotL_8c.html#a058d843996bb36b73784b80aae05f04b" title=" " alt="" coords="403,1525,488,1552"/>
<area shape="rect"  href="ilu__dpanel__dfs_8c.html#a0a3d016444b041668956824248d22439" title=" " alt="" coords="391,1576,499,1603"/>
<area shape="rect"  href="dpanel__bmod_8c.html#a192df249a9fc13ad49bf3f2cd79aba65" title=" " alt="" coords="394,917,497,944"/>
<area shape="rect"  href="ilu__dcolumn__dfs_8c.html#af164b7b553eed616e2ed95144698fe7a" title=" " alt="" coords="386,1019,505,1045"/>
<area shape="rect"  href="dcolumn__bmod_8c.html#a5ca322682f98f276feb3c50b31ca56b8" title=" " alt="" coords="389,1120,502,1147"/>
<area shape="rect"  href="ilu__dcopy__to__ucol_8c.html#ab1802613180b46ffdb7b058a42c38716" title=" " alt="" coords="381,867,510,893"/>
</map>
</div>

</div>
</div>
</div><!-- contents -->
<!-- start footer part -->
<hr class="footer"/><address class="footer"><small>
Generated by &#160;<a href="http://www.doxygen.org/index.html">
<img class="footer" src="doxygen.png" alt="doxygen"/>
</a> 1.8.15
</small></address>
</body>
</html>