This file is indexed.

/usr/share/doc/octave/octave.html/Specialized-Solvers.html is in octave-doc 4.2.2-1ubuntu1.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  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
<!DOCTYPE html PUBLIC "-//W3C//DTD HTML 4.01 Transitional//EN" "http://www.w3.org/TR/html4/loose.dtd">
<html>
<!-- Created by GNU Texinfo 6.5, http://www.gnu.org/software/texinfo/ -->
<head>
<meta http-equiv="Content-Type" content="text/html; charset=utf-8">
<title>Specialized Solvers (GNU Octave)</title>

<meta name="description" content="Specialized Solvers (GNU Octave)">
<meta name="keywords" content="Specialized Solvers (GNU Octave)">
<meta name="resource-type" content="document">
<meta name="distribution" content="global">
<meta name="Generator" content="makeinfo">
<link href="index.html#Top" rel="start" title="Top">
<link href="Concept-Index.html#Concept-Index" rel="index" title="Concept Index">
<link href="index.html#SEC_Contents" rel="contents" title="Table of Contents">
<link href="Linear-Algebra.html#Linear-Algebra" rel="up" title="Linear Algebra">
<link href="Vectorization-and-Faster-Code-Execution.html#Vectorization-and-Faster-Code-Execution" rel="next" title="Vectorization and Faster Code Execution">
<link href="Functions-of-a-Matrix.html#Functions-of-a-Matrix" rel="prev" title="Functions of a Matrix">
<style type="text/css">
<!--
a.summary-letter {text-decoration: none}
blockquote.indentedblock {margin-right: 0em}
blockquote.smallindentedblock {margin-right: 0em; font-size: smaller}
blockquote.smallquotation {font-size: smaller}
div.display {margin-left: 3.2em}
div.example {margin-left: 3.2em}
div.lisp {margin-left: 3.2em}
div.smalldisplay {margin-left: 3.2em}
div.smallexample {margin-left: 3.2em}
div.smalllisp {margin-left: 3.2em}
kbd {font-style: oblique}
pre.display {font-family: inherit}
pre.format {font-family: inherit}
pre.menu-comment {font-family: serif}
pre.menu-preformatted {font-family: serif}
pre.smalldisplay {font-family: inherit; font-size: smaller}
pre.smallexample {font-size: smaller}
pre.smallformat {font-family: inherit; font-size: smaller}
pre.smalllisp {font-size: smaller}
span.nolinebreak {white-space: nowrap}
span.roman {font-family: initial; font-weight: normal}
span.sansserif {font-family: sans-serif; font-weight: normal}
ul.no-bullet {list-style: none}
-->
</style>
<link rel="stylesheet" type="text/css" href="octave.css">


</head>

<body lang="en">
<a name="Specialized-Solvers"></a>
<div class="header">
<p>
Previous: <a href="Functions-of-a-Matrix.html#Functions-of-a-Matrix" accesskey="p" rel="prev">Functions of a Matrix</a>, Up: <a href="Linear-Algebra.html#Linear-Algebra" accesskey="u" rel="up">Linear Algebra</a> &nbsp; [<a href="index.html#SEC_Contents" title="Table of contents" rel="contents">Contents</a>][<a href="Concept-Index.html#Concept-Index" title="Index" rel="index">Index</a>]</p>
</div>
<hr>
<a name="Specialized-Solvers-1"></a>
<h3 class="section">18.5 Specialized Solvers</h3>
<a name="index-matrix_002c-specialized-solvers"></a>

<a name="XREFbicg"></a><dl>
<dt><a name="index-bicg"></a>: <em><var>x</var> =</em> <strong>bicg</strong> <em>(<var>A</var>, <var>b</var>, <var>rtol</var>, <var>maxit</var>, <var>M1</var>, <var>M2</var>, <var>x0</var>)</em></dt>
<dt><a name="index-bicg-1"></a>: <em><var>x</var> =</em> <strong>bicg</strong> <em>(<var>A</var>, <var>b</var>, <var>rtol</var>, <var>maxit</var>, <var>P</var>)</em></dt>
<dt><a name="index-bicg-2"></a>: <em>[<var>x</var>, <var>flag</var>, <var>relres</var>, <var>iter</var>, <var>resvec</var>] =</em> <strong>bicg</strong> <em>(<var>A</var>, <var>b</var>, &hellip;)</em></dt>
<dd><p>Solve <code>A x = b</code> using the Bi-conjugate gradient iterative method.
</p>
<ul class="no-bullet">
<li>- <var>rtol</var> is the relative tolerance, if not given or set to [] the
default value 1e-6 is used.

</li><li>- <var>maxit</var> the maximum number of outer iterations, if not given or
set to [] the default value <code>min (20, numel (b))</code> is used.

</li><li>- <var>x0</var> the initial guess, if not given or set to [] the default
value <code>zeros (size (b))</code> is used.
</li></ul>

<p><var>A</var> can be passed as a matrix or as a function handle or inline function
<code>f</code> such that <code>f(x, &quot;notransp&quot;) = A*x</code> and
<code>f(x, &quot;transp&quot;) = A'*x</code>.
</p>
<p>The preconditioner <var>P</var> is given as <code>P = M1 * M2</code>.  Both <var>M1</var>
and <var>M2</var> can be passed as a matrix or as a function handle or inline
function <code>g</code> such that <code>g(x, &quot;notransp&quot;) = M1 \ x</code> or
<code>g(x, &quot;notransp&quot;) = M2 \ x</code> and <code>g(x, &quot;transp&quot;) = M1' \ x</code> or
<code>g(x, &quot;transp&quot;) = M2' \ x</code>.
</p>
<p>If called with more than one output parameter
</p>
<ul class="no-bullet">
<li>- <var>flag</var> indicates the exit status:

<ul class="no-bullet">
<li>- 0: iteration converged to the within the chosen tolerance

</li><li>- 1: the maximum number of iterations was reached before convergence

</li><li>- 3: the algorithm reached stagnation
</li></ul>

<p>(the value 2 is unused but skipped for compatibility).
</p>
</li><li>- <var>relres</var> is the final value of the relative residual.

</li><li>- <var>iter</var> is the number of iterations performed.

</li><li>- <var>resvec</var> is a vector containing the relative residual at each
iteration.
</li></ul>


<p><strong>See also:</strong> <a href="#XREFbicgstab">bicgstab</a>, <a href="#XREFcgs">cgs</a>, <a href="#XREFgmres">gmres</a>, <a href="Iterative-Techniques.html#XREFpcg">pcg</a>, <a href="#XREFqmr">qmr</a>.
</p>
</dd></dl>


<a name="XREFbicgstab"></a><dl>
<dt><a name="index-bicgstab"></a>: <em><var>x</var> =</em> <strong>bicgstab</strong> <em>(<var>A</var>, <var>b</var>, <var>rtol</var>, <var>maxit</var>, <var>M1</var>, <var>M2</var>, <var>x0</var>)</em></dt>
<dt><a name="index-bicgstab-1"></a>: <em><var>x</var> =</em> <strong>bicgstab</strong> <em>(<var>A</var>, <var>b</var>, <var>rtol</var>, <var>maxit</var>, <var>P</var>)</em></dt>
<dt><a name="index-bicgstab-2"></a>: <em>[<var>x</var>, <var>flag</var>, <var>relres</var>, <var>iter</var>, <var>resvec</var>] =</em> <strong>bicgstab</strong> <em>(<var>A</var>, <var>b</var>, &hellip;)</em></dt>
<dd><p>Solve <code>A x = b</code> using the stabilizied Bi-conjugate gradient iterative
method.
</p>
<ul class="no-bullet">
<li>- <var>rtol</var> is the relative tolerance, if not given or set to [] the
default value 1e-6 is used.

</li><li>- <var>maxit</var> the maximum number of outer iterations, if not given or
set to [] the default value <code>min (20, numel (b))</code> is used.

</li><li>- <var>x0</var> the initial guess, if not given or set to [] the default
value <code>zeros (size (b))</code> is used.
</li></ul>

<p><var>A</var> can be passed as a matrix or as a function handle or inline
function <code>f</code> such that <code>f(x) = A*x</code>.
</p>
<p>The preconditioner <var>P</var> is given as <code>P = M1 * M2</code>.  Both <var>M1</var>
and <var>M2</var> can be passed as a matrix or as a function handle or inline
function <code>g</code> such that <code>g(x) = M1 \ x</code> or <code>g(x) = M2 \ x</code>.
</p>
<p>If called with more than one output parameter
</p>
<ul class="no-bullet">
<li>- <var>flag</var> indicates the exit status:

<ul class="no-bullet">
<li>- 0: iteration converged to the within the chosen tolerance

</li><li>- 1: the maximum number of iterations was reached before convergence

</li><li>- 3: the algorithm reached stagnation
</li></ul>

<p>(the value 2 is unused but skipped for compatibility).
</p>
</li><li>- <var>relres</var> is the final value of the relative residual.

</li><li>- <var>iter</var> is the number of iterations performed.

</li><li>- <var>resvec</var> is a vector containing the relative residual at each
iteration.
</li></ul>


<p><strong>See also:</strong> <a href="#XREFbicg">bicg</a>, <a href="#XREFcgs">cgs</a>, <a href="#XREFgmres">gmres</a>, <a href="Iterative-Techniques.html#XREFpcg">pcg</a>, <a href="#XREFqmr">qmr</a>.
</p>
</dd></dl>


<a name="XREFcgs"></a><dl>
<dt><a name="index-cgs"></a>: <em><var>x</var> =</em> <strong>cgs</strong> <em>(<var>A</var>, <var>b</var>, <var>rtol</var>, <var>maxit</var>, <var>M1</var>, <var>M2</var>, <var>x0</var>)</em></dt>
<dt><a name="index-cgs-1"></a>: <em><var>x</var> =</em> <strong>cgs</strong> <em>(<var>A</var>, <var>b</var>, <var>rtol</var>, <var>maxit</var>, <var>P</var>)</em></dt>
<dt><a name="index-cgs-2"></a>: <em>[<var>x</var>, <var>flag</var>, <var>relres</var>, <var>iter</var>, <var>resvec</var>] =</em> <strong>cgs</strong> <em>(<var>A</var>, <var>b</var>, &hellip;)</em></dt>
<dd><p>Solve <code>A x = b</code>, where <var>A</var> is a square matrix, using the
Conjugate Gradients Squared method.
</p>
<ul class="no-bullet">
<li>- <var>rtol</var> is the relative tolerance, if not given or set to [] the
default value 1e-6 is used.

</li><li>- <var>maxit</var> the maximum number of outer iterations, if not given or
set to [] the default value <code>min (20, numel (b))</code> is used.

</li><li>- <var>x0</var> the initial guess, if not given or set to [] the default
value <code>zeros (size (b))</code> is used.
</li></ul>

<p><var>A</var> can be passed as a matrix or as a function handle or inline
function <code>f</code> such that <code>f(x) = A*x</code>.
</p>
<p>The preconditioner <var>P</var> is given as <code>P = M1 * M2</code>.  Both <var>M1</var>
and <var>M2</var> can be passed as a matrix or as a function handle or inline
function <code>g</code> such that <code>g(x) = M1 \ x</code> or <code>g(x) = M2 \ x</code>.
</p>
<p>If called with more than one output parameter
</p>
<ul class="no-bullet">
<li>- <var>flag</var> indicates the exit status:

<ul class="no-bullet">
<li>- 0: iteration converged to the within the chosen tolerance

</li><li>- 1: the maximum number of iterations was reached before convergence

</li><li>- 3: the algorithm reached stagnation
</li></ul>

<p>(the value 2 is unused but skipped for compatibility).
</p>
</li><li>- <var>relres</var> is the final value of the relative residual.

</li><li>- <var>iter</var> is the number of iterations performed.

</li><li>- <var>resvec</var> is a vector containing the relative residual at
each iteration.
</li></ul>


<p><strong>See also:</strong> <a href="Iterative-Techniques.html#XREFpcg">pcg</a>, <a href="#XREFbicgstab">bicgstab</a>, <a href="#XREFbicg">bicg</a>, <a href="#XREFgmres">gmres</a>, <a href="#XREFqmr">qmr</a>.
</p></dd></dl>


<a name="XREFgmres"></a><dl>
<dt><a name="index-gmres"></a>: <em><var>x</var> =</em> <strong>gmres</strong> <em>(<var>A</var>, <var>b</var>, <var>m</var>, <var>rtol</var>, <var>maxit</var>, <var>M1</var>, <var>M2</var>, <var>x0</var>)</em></dt>
<dt><a name="index-gmres-1"></a>: <em><var>x</var> =</em> <strong>gmres</strong> <em>(<var>A</var>, <var>b</var>, <var>m</var>, <var>rtol</var>, <var>maxit</var>, <var>P</var>)</em></dt>
<dt><a name="index-gmres-2"></a>: <em>[<var>x</var>, <var>flag</var>, <var>relres</var>, <var>iter</var>, <var>resvec</var>] =</em> <strong>gmres</strong> <em>(&hellip;)</em></dt>
<dd><p>Solve <code>A x = b</code> using the Preconditioned GMRES iterative method with
restart, a.k.a. PGMRES(m).
</p>
<ul class="no-bullet">
<li>- <var>rtol</var> is the relative tolerance,
if not given or set to [] the default value 1e-6 is used.

</li><li>- <var>maxit</var> is the maximum number of outer iterations, if not given or
set to [] the default value <code>min (10, numel (b) / restart)</code> is used.

</li><li>- <var>x0</var> is the initial guess,
if not given or set to [] the default value <code>zeros (size (b))</code> is used.

</li><li>- <var>m</var> is the restart parameter,
if not given or set to [] the default value <code>numel (b)</code> is used.
</li></ul>

<p>Argument <var>A</var> can be passed as a matrix, function handle, or inline
function <code>f</code> such that <code>f(x) = A*x</code>.
</p>
<p>The preconditioner <var>P</var> is given as <code>P = M1 * M2</code>.  Both <var>M1</var>
and <var>M2</var> can be passed as a matrix, function handle, or inline function
<code>g</code> such that <code>g(x) = M1\x</code> or <code>g(x) = M2\x</code>.
</p>
<p>Besides the vector <var>x</var>, additional outputs are:
</p>
<ul class="no-bullet">
<li>- <var>flag</var> indicates the exit status:

<dl compact="compact">
<dt>0 : iteration converged to within the specified tolerance</dt>
<dt>1 : maximum number of iterations exceeded</dt>
<dt>2 : unused, but skipped for compatibility</dt>
<dt>3 : algorithm reached stagnation (no change between iterations)</dt>
</dl>

</li><li>- <var>relres</var> is the final value of the relative residual.

</li><li>- <var>iter</var> is a vector containing the number of outer iterations and
total iterations performed.

</li><li>- <var>resvec</var> is a vector containing the relative residual at each
iteration.
</li></ul>


<p><strong>See also:</strong> <a href="#XREFbicg">bicg</a>, <a href="#XREFbicgstab">bicgstab</a>, <a href="#XREFcgs">cgs</a>, <a href="Iterative-Techniques.html#XREFpcg">pcg</a>, <a href="Iterative-Techniques.html#XREFpcr">pcr</a>, <a href="#XREFqmr">qmr</a>.
</p></dd></dl>


<a name="XREFqmr"></a><dl>
<dt><a name="index-qmr"></a>: <em><var>x</var> =</em> <strong>qmr</strong> <em>(<var>A</var>, <var>b</var>, <var>rtol</var>, <var>maxit</var>, <var>M1</var>, <var>M2</var>, <var>x0</var>)</em></dt>
<dt><a name="index-qmr-1"></a>: <em><var>x</var> =</em> <strong>qmr</strong> <em>(<var>A</var>, <var>b</var>, <var>rtol</var>, <var>maxit</var>, <var>P</var>)</em></dt>
<dt><a name="index-qmr-2"></a>: <em>[<var>x</var>, <var>flag</var>, <var>relres</var>, <var>iter</var>, <var>resvec</var>] =</em> <strong>qmr</strong> <em>(<var>A</var>, <var>b</var>, &hellip;)</em></dt>
<dd><p>Solve <code>A x = b</code> using the Quasi-Minimal Residual iterative method
(without look-ahead).
</p>
<ul class="no-bullet">
<li>- <var>rtol</var> is the relative tolerance, if not given or set to [] the
default value 1e-6 is used.

</li><li>- <var>maxit</var> the maximum number of outer iterations, if not given or
set to [] the default value <code>min (20, numel (b))</code> is used.

</li><li>- <var>x0</var> the initial guess, if not given or set to [] the default
value <code>zeros (size (b))</code> is used.
</li></ul>

<p><var>A</var> can be passed as a matrix or as a function handle or inline
function <code>f</code> such that <code>f(x, &quot;notransp&quot;) = A*x</code> and
<code>f(x, &quot;transp&quot;) = A'*x</code>.
</p>
<p>The preconditioner <var>P</var> is given as <code>P = M1 * M2</code>.  Both <var>M1</var>
and <var>M2</var> can be passed as a matrix or as a function handle or inline
function <code>g</code> such that <code>g(x, &quot;notransp&quot;) = M1 \ x</code> or
<code>g(x, &quot;notransp&quot;) = M2 \ x</code> and <code>g(x, &quot;transp&quot;) = M1' \ x</code> or
<code>g(x, &quot;transp&quot;) = M2' \ x</code>.
</p>
<p>If called with more than one output parameter
</p>
<ul class="no-bullet">
<li>- <var>flag</var> indicates the exit status:

<ul class="no-bullet">
<li>- 0: iteration converged to the within the chosen tolerance

</li><li>- 1: the maximum number of iterations was reached before convergence

</li><li>- 3: the algorithm reached stagnation
</li></ul>

<p>(the value 2 is unused but skipped for compatibility).
</p>
</li><li>- <var>relres</var> is the final value of the relative residual.

</li><li>- <var>iter</var> is the number of iterations performed.

</li><li>- <var>resvec</var> is a vector containing the residual norms at each
      iteration.
</li></ul>

<p>References:
</p>
<ol>
<li> R. Freund and N. Nachtigal, <cite>QMR: a quasi-minimal residual
method for non-Hermitian linear systems</cite>, Numerische Mathematik,
1991, 60, pp. 315-339.

</li><li> R. Barrett, M. Berry, T. Chan, J. Demmel, J. Donato, J. Dongarra,
V. Eijkhour, R. Pozo, C. Romine, and H. van der Vorst,
<cite>Templates for the solution of linear systems: Building blocks
for iterative methods</cite>, SIAM, 2nd ed., 1994.
</li></ol>


<p><strong>See also:</strong> <a href="#XREFbicg">bicg</a>, <a href="#XREFbicgstab">bicgstab</a>, <a href="#XREFcgs">cgs</a>, <a href="#XREFgmres">gmres</a>, <a href="Iterative-Techniques.html#XREFpcg">pcg</a>.
</p></dd></dl>



<hr>
<div class="header">
<p>
Previous: <a href="Functions-of-a-Matrix.html#Functions-of-a-Matrix" accesskey="p" rel="prev">Functions of a Matrix</a>, Up: <a href="Linear-Algebra.html#Linear-Algebra" accesskey="u" rel="up">Linear Algebra</a> &nbsp; [<a href="index.html#SEC_Contents" title="Table of contents" rel="contents">Contents</a>][<a href="Concept-Index.html#Concept-Index" title="Index" rel="index">Index</a>]</p>
</div>



</body>
</html>