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 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876
|
<!DOCTYPE html>
<html>
<!-- Created by GNU Texinfo 7.1.1, https://www.gnu.org/software/texinfo/ -->
<head>
<meta http-equiv="Content-Type" content="text/html; charset=utf-8">
<title>Basic Matrix Functions (GNU Octave (version 10.3.0))</title>
<meta name="description" content="Basic Matrix Functions (GNU Octave (version 10.3.0))">
<meta name="keywords" content="Basic Matrix Functions (GNU Octave (version 10.3.0))">
<meta name="resource-type" content="document">
<meta name="distribution" content="global">
<meta name="Generator" content="makeinfo">
<meta name="viewport" content="width=device-width,initial-scale=1">
<link href="index.html" rel="start" title="Top">
<link href="Concept-Index.html" rel="index" title="Concept Index">
<link href="index.html#SEC_Contents" rel="contents" title="Table of Contents">
<link href="Linear-Algebra.html" rel="up" title="Linear Algebra">
<link href="Matrix-Factorizations.html" rel="next" title="Matrix Factorizations">
<link href="Techniques-Used-for-Linear-Algebra.html" rel="prev" title="Techniques Used for Linear Algebra">
<style type="text/css">
<!--
a.copiable-link {visibility: hidden; text-decoration: none; line-height: 0em}
div.example {margin-left: 3.2em}
span:hover a.copiable-link {visibility: visible}
strong.def-name {font-family: monospace; font-weight: bold; font-size: larger}
ul.mark-bullet {list-style-type: disc}
-->
</style>
<link rel="stylesheet" type="text/css" href="octave.css">
</head>
<body lang="en">
<div class="section-level-extent" id="Basic-Matrix-Functions">
<div class="nav-panel">
<p>
Next: <a href="Matrix-Factorizations.html" accesskey="n" rel="next">Matrix Factorizations</a>, Previous: <a href="Techniques-Used-for-Linear-Algebra.html" accesskey="p" rel="prev">Techniques Used for Linear Algebra</a>, Up: <a href="Linear-Algebra.html" accesskey="u" rel="up">Linear Algebra</a> [<a href="index.html#SEC_Contents" title="Table of contents" rel="contents">Contents</a>][<a href="Concept-Index.html" title="Index" rel="index">Index</a>]</p>
</div>
<hr>
<h3 class="section" id="Basic-Matrix-Functions-1"><span>18.2 Basic Matrix Functions<a class="copiable-link" href="#Basic-Matrix-Functions-1"> ¶</a></span></h3>
<a class="index-entry-id" id="index-matrix-functions_002c-basic"></a>
<a class="anchor" id="XREFbalance"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-balance"><span><code class="def-type"><var class="var">AA</var> =</code> <strong class="def-name">balance</strong> <code class="def-code-arguments">(<var class="var">A</var>)</code><a class="copiable-link" href="#index-balance"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-balance-1"><span><code class="def-type"><var class="var">AA</var> =</code> <strong class="def-name">balance</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">opt</var>)</code><a class="copiable-link" href="#index-balance-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-balance-2"><span><code class="def-type">[<var class="var">DD</var>, <var class="var">AA</var>] =</code> <strong class="def-name">balance</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">opt</var>)</code><a class="copiable-link" href="#index-balance-2"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-balance-3"><span><code class="def-type">[<var class="var">D</var>, <var class="var">P</var>, <var class="var">AA</var>] =</code> <strong class="def-name">balance</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">opt</var>)</code><a class="copiable-link" href="#index-balance-3"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-balance-4"><span><code class="def-type">[<var class="var">CC</var>, <var class="var">DD</var>, <var class="var">AA</var>, <var class="var">BB</var>] =</code> <strong class="def-name">balance</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">B</var>, <var class="var">opt</var>)</code><a class="copiable-link" href="#index-balance-4"> ¶</a></span></dt>
<dd>
<p>Balance the matrix <var class="var">A</var> to reduce numerical errors in future
calculations.
</p>
<p>Compute <code class="code"><var class="var">AA</var> = <var class="var">DD</var> \ <var class="var">A</var> * <var class="var">DD</var></code> in which <var class="var">AA</var>
is a matrix whose row and column norms are roughly equal in magnitude, and
<code class="code"><var class="var">DD</var> = <var class="var">P</var> * <var class="var">D</var></code>, in which <var class="var">P</var> is a permutation
matrix and <var class="var">D</var> is a diagonal matrix of powers of two. This allows the
equilibration to be computed without round-off. Results of eigenvalue
calculation are typically improved by balancing first.
</p>
<p>If two output values are requested, <code class="code">balance</code> returns
the diagonal <var class="var">D</var> and the permutation <var class="var">P</var> separately as vectors.
In this case, <code class="code"><var class="var">DD</var> = eye(n)(:,<var class="var">P</var>) * diag (<var class="var">D</var>)</code>, where
<em class="math">n</em> is the matrix size.
</p>
<p>If four output values are requested, compute <code class="code"><var class="var">AA</var> =
<var class="var">CC</var>*<var class="var">A</var>*<var class="var">DD</var></code> and <code class="code"><var class="var">BB</var> = <var class="var">CC</var>*<var class="var">B</var>*<var class="var">DD</var></code>,
in which <var class="var">AA</var> and <var class="var">BB</var> have nonzero elements of approximately the
same magnitude and <var class="var">CC</var> and <var class="var">DD</var> are permuted diagonal matrices as
in <var class="var">DD</var> for the algebraic eigenvalue problem.
</p>
<p>The eigenvalue balancing option <var class="var">opt</var> may be one of:
</p>
<dl class="table">
<dt><code class="code">"noperm"</code>, <code class="code">"S"</code></dt>
<dd><p>Scale only; do not permute.
</p>
</dd>
<dt><code class="code">"noscal"</code>, <code class="code">"P"</code></dt>
<dd><p>Permute only; do not scale.
</p></dd>
</dl>
<p>Algebraic eigenvalue balancing uses standard <small class="sc">LAPACK</small> routines.
</p>
<p>Generalized eigenvalue problem balancing uses Ward’s algorithm
(SIAM Journal on Scientific and Statistical Computing, 1981).
</p></dd></dl>
<a class="anchor" id="XREFbandwidth"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-bandwidth"><span><code class="def-type"><var class="var">bw</var> =</code> <strong class="def-name">bandwidth</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">type</var>)</code><a class="copiable-link" href="#index-bandwidth"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-bandwidth-1"><span><code class="def-type">[<var class="var">lower</var>, <var class="var">upper</var>] =</code> <strong class="def-name">bandwidth</strong> <code class="def-code-arguments">(<var class="var">A</var>)</code><a class="copiable-link" href="#index-bandwidth-1"> ¶</a></span></dt>
<dd><p>Compute the bandwidth of <var class="var">A</var>.
</p>
<p>The <var class="var">type</var> argument is the string <code class="code">"lower"</code> for the lower
bandwidth and <code class="code">"upper"</code> for the upper bandwidth. If no <var class="var">type</var> is
specified return both the lower and upper bandwidth of <var class="var">A</var>.
</p>
<p>The lower/upper bandwidth of a matrix is the number of
subdiagonals/superdiagonals with nonzero entries.
</p>
<p><strong class="strong">See also:</strong> <a class="ref" href="Predicates-for-Numeric-Objects.html#XREFisbanded">isbanded</a>, <a class="ref" href="Predicates-for-Numeric-Objects.html#XREFisdiag">isdiag</a>, <a class="ref" href="Predicates-for-Numeric-Objects.html#XREFistril">istril</a>, <a class="ref" href="Predicates-for-Numeric-Objects.html#XREFistriu">istriu</a>.
</p></dd></dl>
<a class="anchor" id="XREFcond"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-cond"><span><code class="def-type"><var class="var">c</var> =</code> <strong class="def-name">cond</strong> <code class="def-code-arguments">(<var class="var">A</var>)</code><a class="copiable-link" href="#index-cond"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-cond-1"><span><code class="def-type"><var class="var">c</var> =</code> <strong class="def-name">cond</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">p</var>)</code><a class="copiable-link" href="#index-cond-1"> ¶</a></span></dt>
<dd><p>Compute the <var class="var">p</var>-norm condition number of a matrix with respect to
inversion.
</p>
<p><code class="code">cond (<var class="var">A</var>)</code> is defined as
<code class="code">norm (<var class="var">A</var>, <var class="var">p</var>) * norm (inv (<var class="var">A</var>), <var class="var">p</var>)</code>.
</p>
<p>By default, <code class="code"><var class="var">p</var> = 2</code> is used which implies a (relatively slow)
singular value decomposition. Other possible selections are
<code class="code"><var class="var">p</var> = 1, Inf, "fro"</code> which are generally faster. For a full
discussion of possible <var class="var">p</var> values, see <a class="pxref" href="#XREFnorm"><code class="code">norm</code></a>.
</p>
<p>The condition number of a matrix quantifies the sensitivity of the matrix
inversion operation when small changes are made to matrix elements. Ideally
the condition number will be close to 1. When the number is large this
indicates small changes (such as underflow or round-off error) will produce
large changes in the resulting output. In such cases the solution results
from numerical computing are not likely to be accurate.
</p>
<p><strong class="strong">See also:</strong> <a class="ref" href="Sparse-Linear-Algebra.html#XREFcondest">condest</a>, <a class="ref" href="#XREFrcond">rcond</a>, <a class="ref" href="#XREFcondeig">condeig</a>, <a class="ref" href="#XREFnorm">norm</a>, <a class="ref" href="Matrix-Factorizations.html#XREFsvd">svd</a>.
</p></dd></dl>
<a class="anchor" id="XREFcondeig"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-condeig"><span><code class="def-type"><var class="var">c</var> =</code> <strong class="def-name">condeig</strong> <code class="def-code-arguments">(<var class="var">a</var>)</code><a class="copiable-link" href="#index-condeig"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-condeig-1"><span><code class="def-type">[<var class="var">v</var>, <var class="var">lambda</var>, <var class="var">c</var>] =</code> <strong class="def-name">condeig</strong> <code class="def-code-arguments">(<var class="var">a</var>)</code><a class="copiable-link" href="#index-condeig-1"> ¶</a></span></dt>
<dd><p>Compute condition numbers of a matrix with respect to eigenvalues.
</p>
<p>The condition numbers are the reciprocals of the cosines of the angles
between the left and right eigenvectors; Large values indicate that the
matrix has multiple distinct eigenvalues.
</p>
<p>The input <var class="var">a</var> must be a square numeric matrix.
</p>
<p>The outputs are:
</p>
<ul class="itemize mark-bullet">
<li><var class="var">c</var> is a vector of condition numbers for the eigenvalues of
<var class="var">a</var>.
</li><li><var class="var">v</var> is the matrix of right eigenvectors of <var class="var">a</var>. The result is
equivalent to calling <code class="code">[<var class="var">v</var>, <var class="var">lambda</var>] = eig (<var class="var">a</var>)</code>.
</li><li><var class="var">lambda</var> is the diagonal matrix of eigenvalues of <var class="var">a</var>. The
result is equivalent to calling
<code class="code">[<var class="var">v</var>, <var class="var">lambda</var>] = eig (<var class="var">a</var>)</code>.
</li></ul>
<p>Example
</p>
<div class="example">
<div class="group"><pre class="example-preformatted">a = [1, 2; 3, 4];
c = condeig (a)
⇒ c =
1.0150
1.0150
</pre></div></div>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFeig">eig</a>, <a class="ref" href="#XREFcond">cond</a>, <a class="ref" href="#XREFbalance">balance</a>.
</p></dd></dl>
<a class="anchor" id="XREFdet"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-det"><span><code class="def-type"><var class="var">d</var> =</code> <strong class="def-name">det</strong> <code class="def-code-arguments">(<var class="var">A</var>)</code><a class="copiable-link" href="#index-det"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-det-1"><span><code class="def-type">[<var class="var">d</var>, <var class="var">rcond</var>] =</code> <strong class="def-name">det</strong> <code class="def-code-arguments">(<var class="var">A</var>)</code><a class="copiable-link" href="#index-det-1"> ¶</a></span></dt>
<dd><p>Compute the determinant of <var class="var">A</var>.
</p>
<p>Return an estimate of the reciprocal condition number if requested.
</p>
<p>Programming Notes: Routines from <small class="sc">LAPACK</small> are used for full matrices and
code from <small class="sc">UMFPACK</small> is used for sparse matrices.
</p>
<p>The determinant should not be used to check a matrix for singularity.
For that, use any of the condition number functions: <code class="code">cond</code>,
<code class="code">condest</code>, <code class="code">rcond</code>.
</p>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFcond">cond</a>, <a class="ref" href="Sparse-Linear-Algebra.html#XREFcondest">condest</a>, <a class="ref" href="#XREFrcond">rcond</a>.
</p></dd></dl>
<a class="anchor" id="XREFeig"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-eig"><span><code class="def-type"><var class="var">lambda</var> =</code> <strong class="def-name">eig</strong> <code class="def-code-arguments">(<var class="var">A</var>)</code><a class="copiable-link" href="#index-eig"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-eig-1"><span><code class="def-type"><var class="var">lambda</var> =</code> <strong class="def-name">eig</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">B</var>)</code><a class="copiable-link" href="#index-eig-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-eig-2"><span><code class="def-type">[<var class="var">V</var>, <var class="var">lambda</var>] =</code> <strong class="def-name">eig</strong> <code class="def-code-arguments">(<var class="var">A</var>)</code><a class="copiable-link" href="#index-eig-2"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-eig-3"><span><code class="def-type">[<var class="var">V</var>, <var class="var">lambda</var>] =</code> <strong class="def-name">eig</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">B</var>)</code><a class="copiable-link" href="#index-eig-3"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-eig-4"><span><code class="def-type">[<var class="var">V</var>, <var class="var">lambda</var>, <var class="var">W</var>] =</code> <strong class="def-name">eig</strong> <code class="def-code-arguments">(<var class="var">A</var>)</code><a class="copiable-link" href="#index-eig-4"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-eig-5"><span><code class="def-type">[<var class="var">V</var>, <var class="var">lambda</var>, <var class="var">W</var>] =</code> <strong class="def-name">eig</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">B</var>)</code><a class="copiable-link" href="#index-eig-5"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-eig-6"><span><code class="def-type">[…] =</code> <strong class="def-name">eig</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">balanceOption</var>)</code><a class="copiable-link" href="#index-eig-6"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-eig-7"><span><code class="def-type">[…] =</code> <strong class="def-name">eig</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">B</var>, <var class="var">algorithm</var>)</code><a class="copiable-link" href="#index-eig-7"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-eig-8"><span><code class="def-type">[…] =</code> <strong class="def-name">eig</strong> <code class="def-code-arguments">(…, <var class="var">eigvalOption</var>)</code><a class="copiable-link" href="#index-eig-8"> ¶</a></span></dt>
<dd><p>Compute the eigenvalues (<var class="var">lambda</var>) and optionally the right eigenvectors
(<var class="var">V</var>) and the left eigenvectors (<var class="var">W</var>) of a matrix or pair of matrices.
</p>
<p>The flag <var class="var">balanceOption</var> can be one of:
</p>
<dl class="table">
<dt><code class="code">"balance"</code> (default)</dt>
<dd><p>Preliminary balancing is on.
</p>
</dd>
<dt><code class="code">"nobalance"</code></dt>
<dd><p>Disables preliminary balancing.
</p></dd>
</dl>
<p>The flag <var class="var">eigvalOption</var> can be one of:
</p>
<dl class="table">
<dt><code class="code">"matrix"</code></dt>
<dd><p>Return the eigenvalues in a diagonal matrix. (default if 2 or 3 outputs
are requested)
</p>
</dd>
<dt><code class="code">"vector"</code></dt>
<dd><p>Return the eigenvalues in a column vector. (default if only 1 output is
requested, e.g., <code class="code"><var class="var">lambda</var> = eig (<var class="var">A</var>)</code>)
</p></dd>
</dl>
<p>The flag <var class="var">algorithm</var> can be one of:
</p>
<dl class="table">
<dt><code class="code">"chol"</code></dt>
<dd><p>Use the Cholesky factorization of B. (default if <var class="var">A</var> is symmetric
(Hermitian) and <var class="var">B</var> is symmetric (Hermitian) positive definite)
</p>
</dd>
<dt><code class="code">"qz"</code></dt>
<dd><p>Use the QZ algorithm. (used whenever <var class="var">A</var> or <var class="var">B</var> are not symmetric)
</p></dd>
</dl>
<table class="multitable">
<thead><tr><th width="44%">A and B</th><th width="14%">no flag</th><th width="14%">chol</th><th width="10%">qz</th></tr></thead>
<tbody><tr><td width="44%">both are symmetric</td><td width="14%"><code class="code">"chol"</code></td><td width="14%"><code class="code">"chol"</code></td><td width="10%"><code class="code">"qz"</code></td></tr>
<tr><td width="44%">at least one is not symmetric</td><td width="14%"><code class="code">"qz"</code></td><td width="14%"><code class="code">"qz"</code></td><td width="10%"><code class="code">"qz"</code></td></tr>
</tbody>
</table>
<p>The eigenvalues returned by <code class="code">eig</code> are not ordered.
</p>
<p><strong class="strong">See also:</strong> <a class="ref" href="Sparse-Linear-Algebra.html#XREFeigs">eigs</a>, <a class="ref" href="Matrix-Factorizations.html#XREFsvd">svd</a>.
</p></dd></dl>
<a class="anchor" id="XREFgivens"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-givens"><span><code class="def-type"><var class="var">G</var> =</code> <strong class="def-name">givens</strong> <code class="def-code-arguments">(<var class="var">x</var>, <var class="var">y</var>)</code><a class="copiable-link" href="#index-givens"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-givens-1"><span><code class="def-type">[<var class="var">c</var>, <var class="var">s</var>] =</code> <strong class="def-name">givens</strong> <code class="def-code-arguments">(<var class="var">x</var>, <var class="var">y</var>)</code><a class="copiable-link" href="#index-givens-1"> ¶</a></span></dt>
<dd><p>Compute the Givens rotation matrix <var class="var">G</var>.
The Givens matrix is a 2-by-2 orthogonal matrix
</p>
<div class="example">
<div class="group"><pre class="example-preformatted"><var class="var">G</var> = [ <var class="var">c</var> , <var class="var">s</var>
-<var class="var">s</var>', <var class="var">c</var>]
</pre></div></div>
<p>such that
</p>
<div class="example">
<pre class="example-preformatted"><var class="var">G</var> * [<var class="var">x</var>; <var class="var">y</var>] = [*; 0]
</pre></div>
<p>with <var class="var">x</var> and <var class="var">y</var> scalars.
</p>
<p>If two output arguments are requested, return the factors <var class="var">c</var> and <var class="var">s</var>
rather than the Givens rotation matrix.
</p>
<p>For example:
</p>
<div class="example">
<div class="group"><pre class="example-preformatted">givens (1, 1)
⇒ 0.70711 0.70711
-0.70711 0.70711
</pre></div></div>
<p>Note: The Givens matrix represents a counterclockwise rotation of a 2-D
plane and can be used to introduce zeros into a matrix prior to complete
factorization.
</p>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFplanerot">planerot</a>, <a class="ref" href="Matrix-Factorizations.html#XREFqr">qr</a>.
</p></dd></dl>
<a class="anchor" id="XREFgsvd"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-gsvd"><span><code class="def-type"><var class="var">S</var> =</code> <strong class="def-name">gsvd</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">B</var>)</code><a class="copiable-link" href="#index-gsvd"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-gsvd-1"><span><code class="def-type">[<var class="var">U</var>, <var class="var">V</var>, <var class="var">X</var>, <var class="var">C</var>, <var class="var">S</var>] =</code> <strong class="def-name">gsvd</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">B</var>)</code><a class="copiable-link" href="#index-gsvd-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-gsvd-2"><span><code class="def-type">[<var class="var">U</var>, <var class="var">V</var>, <var class="var">X</var>, <var class="var">C</var>, <var class="var">S</var>] =</code> <strong class="def-name">gsvd</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">B</var>, 0)</code><a class="copiable-link" href="#index-gsvd-2"> ¶</a></span></dt>
<dd><p>Compute the generalized singular value decomposition of (<var class="var">A</var>, <var class="var">B</var>).
</p>
<p>The generalized singular value decomposition is defined by the following
relations:
</p>
<div class="example">
<div class="group"><pre class="example-preformatted">A = U*C*X'
B = V*S*X'
C'*C + S'*S = eye (columns (A))
</pre></div></div>
<p>The function <code class="code">gsvd</code> normally returns just the vector of generalized
singular values
<code class="code">sqrt (diag (C'*C) ./ diag (S'*S))</code>.
If asked for five return values, it also computes
U, V, X, and C.
</p>
<p>If the optional third input is present, <code class="code">gsvd</code> constructs the
"economy-sized" decomposition where the number of columns of <var class="var">U</var>, <var class="var">V</var>
and the number of rows of <var class="var">C</var>, <var class="var">S</var> is less than or equal to the number
of columns of <var class="var">A</var>. This option is not yet implemented.
</p>
<p>Programming Note: the code is a wrapper to the corresponding <small class="sc">LAPACK</small> dggsvd
and zggsvd routines. If matrices <var class="var">A</var> and <var class="var">B</var> are <em class="emph">both</em> rank
deficient then <small class="sc">LAPACK</small> will return an incorrect factorization. Programmers
should avoid this combination.
</p>
<p><strong class="strong">See also:</strong> <a class="ref" href="Matrix-Factorizations.html#XREFsvd">svd</a>.
</p></dd></dl>
<a class="anchor" id="XREFplanerot"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-planerot"><span><code class="def-type">[<var class="var">G</var>, <var class="var">y</var>] =</code> <strong class="def-name">planerot</strong> <code class="def-code-arguments">(<var class="var">x</var>)</code><a class="copiable-link" href="#index-planerot"> ¶</a></span></dt>
<dd><p>Compute the Givens rotation matrix for the two-element column vector
<var class="var">x</var>.
The Givens matrix is a 2-by-2 orthogonal matrix
</p>
<div class="example">
<div class="group"><pre class="example-preformatted"><var class="var">G</var> = [ <var class="var">c</var> , <var class="var">s</var>
-<var class="var">s</var>', <var class="var">c</var>]
</pre></div></div>
<p>such that
</p>
<div class="example">
<pre class="example-preformatted"><var class="var">y</var> = <var class="var">G</var> * [<var class="var">x</var>(1); <var class="var">x</var>(2)] ≡ [*; 0]
</pre></div>
<p>Note: The Givens matrix represents a counterclockwise rotation of a 2-D
plane and can be used to introduce zeros into a matrix prior to complete
factorization.
</p>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFgivens">givens</a>, <a class="ref" href="Matrix-Factorizations.html#XREFqr">qr</a>.
</p></dd></dl>
<a class="anchor" id="XREFinv"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-inv"><span><code class="def-type"><var class="var">x</var> =</code> <strong class="def-name">inv</strong> <code class="def-code-arguments">(<var class="var">A</var>)</code><a class="copiable-link" href="#index-inv"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-inv-1"><span><code class="def-type">[<var class="var">x</var>, <var class="var">rcond</var>] =</code> <strong class="def-name">inv</strong> <code class="def-code-arguments">(<var class="var">A</var>)</code><a class="copiable-link" href="#index-inv-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-inverse"><span><code class="def-type">[…] =</code> <strong class="def-name">inverse</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-inverse"> ¶</a></span></dt>
<dd><p>Compute the inverse of the square matrix <var class="var">A</var>.
</p>
<p>Return an estimate of the reciprocal condition number if requested,
otherwise warn of an ill-conditioned matrix if the reciprocal condition
number is small.
</p>
<p>In general it is best to avoid calculating the inverse of a matrix directly.
For example, it is both faster and more accurate to solve systems of
equations (<var class="var">A</var>*<em class="math">x</em> = <em class="math">b</em>) with
<code class="code"><var class="var">y</var> = <var class="var">A</var> \ <em class="math">b</em></code>, rather than
<code class="code"><var class="var">y</var> = inv (<var class="var">A</var>) * <em class="math">b</em></code>.
</p>
<p>If called with a sparse matrix, then in general <var class="var">x</var> will be a full
matrix requiring significantly more storage. Avoid forming the inverse of a
sparse matrix if possible.
</p>
<p>Programming Note: <code class="code">inverse</code> is an alias for <code class="code">inv</code> and can be used
interchangeably.
</p>
<p><strong class="strong">See also:</strong> <a class="ref" href="Arithmetic-Ops.html#XREFldivide">ldivide</a>, <a class="ref" href="Arithmetic-Ops.html#XREFrdivide">rdivide</a>, <a class="ref" href="#XREFpinv">pinv</a>.
</p></dd></dl>
<a class="anchor" id="XREFlinsolve"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-linsolve"><span><code class="def-type"><var class="var">x</var> =</code> <strong class="def-name">linsolve</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">b</var>)</code><a class="copiable-link" href="#index-linsolve"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-linsolve-1"><span><code class="def-type"><var class="var">x</var> =</code> <strong class="def-name">linsolve</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">b</var>, <var class="var">opts</var>)</code><a class="copiable-link" href="#index-linsolve-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-linsolve-2"><span><code class="def-type">[<var class="var">x</var>, <var class="var">R</var>] =</code> <strong class="def-name">linsolve</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-linsolve-2"> ¶</a></span></dt>
<dd><p>Solve the linear system <code class="code">A*x = b</code>.
</p>
<p>With no options, this function is equivalent to the left division operator
(<code class="code">x = A \ b</code>)<!-- /@w --> or the matrix-left-divide function
(<code class="code">x = mldivide (A, b)</code>)<!-- /@w -->.
</p>
<p>Octave ordinarily examines the properties of the matrix <var class="var">A</var> and chooses
a solver that best matches the matrix. By passing a structure <var class="var">opts</var>
to <code class="code">linsolve</code> you can inform Octave directly about the matrix <var class="var">A</var>.
In this case Octave will skip the matrix examination and proceed directly
to solving the linear system.
</p>
<p><strong class="strong">Warning:</strong> If the matrix <var class="var">A</var> does not have the properties listed
in the <var class="var">opts</var> structure then the result will not be accurate AND no
warning will be given. When in doubt, let Octave examine the matrix and
choose the appropriate solver as this step takes little time and the result
is cached so that it is only done once per linear system.
</p>
<p>Possible <var class="var">opts</var> fields (set value to true/false):
</p>
<dl class="table">
<dt>LT</dt>
<dd><p><var class="var">A</var> is lower triangular
</p>
</dd>
<dt>UT</dt>
<dd><p><var class="var">A</var> is upper triangular
</p>
</dd>
<dt>UHESS</dt>
<dd><p><var class="var">A</var> is upper Hessenberg (currently makes no difference)
</p>
</dd>
<dt>SYM</dt>
<dd><p><var class="var">A</var> is symmetric or complex Hermitian (currently makes no difference)
</p>
</dd>
<dt>POSDEF</dt>
<dd><p><var class="var">A</var> is positive definite
</p>
</dd>
<dt>RECT</dt>
<dd><p><var class="var">A</var> is general rectangular (currently makes no difference)
</p>
</dd>
<dt>TRANSA</dt>
<dd><p>Solve <code class="code">A'*x = b</code> if true rather than <code class="code">A*x = b</code>
</p></dd>
</dl>
<p>The optional second output <var class="var">R</var> is the inverse condition number of
<var class="var">A</var> (zero if matrix is singular).
</p>
<p><strong class="strong">See also:</strong> <a class="ref" href="Arithmetic-Ops.html#XREFmldivide">mldivide</a>, <a class="ref" href="#XREFmatrix_005ftype">matrix_type</a>, <a class="ref" href="#XREFrcond">rcond</a>.
</p></dd></dl>
<a class="anchor" id="XREFmatrix_005ftype"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-matrix_005ftype"><span><code class="def-type"><var class="var">type</var> =</code> <strong class="def-name">matrix_type</strong> <code class="def-code-arguments">(<var class="var">A</var>)</code><a class="copiable-link" href="#index-matrix_005ftype"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-matrix_005ftype-1"><span><code class="def-type"><var class="var">type</var> =</code> <strong class="def-name">matrix_type</strong> <code class="def-code-arguments">(<var class="var">A</var>, "nocompute")</code><a class="copiable-link" href="#index-matrix_005ftype-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-matrix_005ftype-2"><span><code class="def-type"><var class="var">A</var> =</code> <strong class="def-name">matrix_type</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">type</var>)</code><a class="copiable-link" href="#index-matrix_005ftype-2"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-matrix_005ftype-3"><span><code class="def-type"><var class="var">A</var> =</code> <strong class="def-name">matrix_type</strong> <code class="def-code-arguments">(<var class="var">A</var>, "upper", <var class="var">perm</var>)</code><a class="copiable-link" href="#index-matrix_005ftype-3"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-matrix_005ftype-4"><span><code class="def-type"><var class="var">A</var> =</code> <strong class="def-name">matrix_type</strong> <code class="def-code-arguments">(<var class="var">A</var>, "lower", <var class="var">perm</var>)</code><a class="copiable-link" href="#index-matrix_005ftype-4"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-matrix_005ftype-5"><span><code class="def-type"><var class="var">A</var> =</code> <strong class="def-name">matrix_type</strong> <code class="def-code-arguments">(<var class="var">A</var>, "banded", <var class="var">nl</var>, <var class="var">nu</var>)</code><a class="copiable-link" href="#index-matrix_005ftype-5"> ¶</a></span></dt>
<dd><p>Identify the matrix type or mark a matrix as a particular type.
</p>
<p>This allows more rapid solutions of linear equations involving <var class="var">A</var> to be
performed.
</p>
<p>Called with a single argument, <code class="code">matrix_type</code> returns the type of the
matrix and caches it for future use.
</p>
<p>Called with more than one argument, <code class="code">matrix_type</code> allows the type of
the matrix to be defined.
</p>
<p>If the option <code class="code">"nocompute"</code> is given, the function will not attempt
to guess the type if it is still unknown. This is useful for debugging
purposes.
</p>
<p>The possible matrix types depend on whether the matrix is full or sparse,
and can be one of the following
</p>
<dl class="table">
<dt><code class="code">"unknown"</code></dt>
<dd><p>Remove any previously cached matrix type, and mark type as unknown.
</p>
</dd>
<dt><code class="code">"full"</code></dt>
<dd><p>Mark the matrix as full.
</p>
</dd>
<dt><code class="code">"positive definite"</code></dt>
<dd><p>Probable full positive definite matrix.
</p>
</dd>
<dt><code class="code">"diagonal"</code></dt>
<dd><p>Diagonal matrix. (Sparse matrices only)
</p>
</dd>
<dt><code class="code">"permuted diagonal"</code></dt>
<dd><p>Permuted Diagonal matrix. The permutation does not need to be specifically
indicated, as the structure of the matrix explicitly gives this. (Sparse
matrices only)
</p>
</dd>
<dt><code class="code">"upper"</code></dt>
<dd><p>Upper triangular. If the optional third argument <var class="var">perm</var> is given, the
matrix is assumed to be a permuted upper triangular with the permutations
defined by the vector <var class="var">perm</var>.
</p>
</dd>
<dt><code class="code">"lower"</code></dt>
<dd><p>Lower triangular. If the optional third argument <var class="var">perm</var> is given, the
matrix is assumed to be a permuted lower triangular with the permutations
defined by the vector <var class="var">perm</var>.
</p>
</dd>
<dt><code class="code">"banded"</code></dt>
<dt><code class="code">"banded positive definite"</code></dt>
<dd><p>Banded matrix with the band size of <var class="var">nl</var> below the diagonal and <var class="var">nu</var>
above it. If <var class="var">nl</var> and <var class="var">nu</var> are 1, then the matrix is tridiagonal
and treated with specialized code. In addition the matrix can be marked as
probably a positive definite. (Sparse matrices only)
</p>
</dd>
<dt><code class="code">"singular"</code></dt>
<dd><p>The matrix is assumed to be singular and will be treated with a minimum norm
solution.
</p>
</dd>
</dl>
<p>Note that the matrix type will be discovered automatically on the first
attempt to solve a linear equation involving <var class="var">A</var>. Therefore
<code class="code">matrix_type</code> is only useful to give Octave hints of the matrix type.
Incorrectly defining the matrix type will result in incorrect results from
solutions of linear equations; it is entirely <strong class="strong">the responsibility of
the user</strong> to correctly identify the matrix type.
</p>
<p>Also, the test for positive definiteness is a low-cost test for a Hermitian
matrix with a real positive diagonal. This does not guarantee that the
matrix is positive definite, but only that it is a probable candidate. When
such a matrix is factorized, a Cholesky factorization is first
attempted, and if that fails the matrix is then treated with an
LU factorization. Once the matrix has been factorized,
<code class="code">matrix_type</code> will return the correct classification of the matrix.
</p></dd></dl>
<a class="anchor" id="XREFnorm"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-norm"><span><code class="def-type"><var class="var">n</var> =</code> <strong class="def-name">norm</strong> <code class="def-code-arguments">(<var class="var">A</var>)</code><a class="copiable-link" href="#index-norm"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-norm-1"><span><code class="def-type"><var class="var">n</var> =</code> <strong class="def-name">norm</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">p</var>)</code><a class="copiable-link" href="#index-norm-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-norm-2"><span><code class="def-type"><var class="var">n</var> =</code> <strong class="def-name">norm</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">p</var>, <var class="var">opt</var>)</code><a class="copiable-link" href="#index-norm-2"> ¶</a></span></dt>
<dd><p>Compute the p-norm of the matrix <var class="var">A</var>.
</p>
<p>If the second argument is not given, <code class="code">p = 2</code><!-- /@w --> is used.
</p>
<p>If <var class="var">A</var> is a matrix (or sparse matrix):
</p>
<dl class="table">
<dt><var class="var">p</var> = <code class="code">1</code></dt>
<dd><p>1-norm, the largest column sum of the absolute values of <var class="var">A</var>.
</p>
</dd>
<dt><var class="var">p</var> = <code class="code">2</code></dt>
<dd><p>Largest singular value of <var class="var">A</var>.
</p>
</dd>
<dt><a id="index-infinity-norm"></a><span><var class="var">p</var> = <code class="code">Inf</code> or <code class="code">"inf"</code><a class="copiable-link" href="#index-infinity-norm"> ¶</a></span></dt>
<dd><p>Infinity norm, the largest row sum of the absolute values of <var class="var">A</var>.
</p>
</dd>
<dt><a id="index-Frobenius-norm"></a><span><var class="var">p</var> = <code class="code">"fro"</code><a class="copiable-link" href="#index-Frobenius-norm"> ¶</a></span></dt>
<dd><p>Frobenius norm of <var class="var">A</var>,
<code class="code">sqrt (sum (diag (<var class="var">A</var>' * <var class="var">A</var>)))</code>.
</p>
</dd>
<dt><a id="index-general-p_002dnorm"></a><span>other <var class="var">p</var>, <code class="code"><var class="var">p</var> > 1</code><a class="copiable-link" href="#index-general-p_002dnorm"> ¶</a></span></dt>
<dd><p>maximum <code class="code">norm (A*x, p)</code> such that <code class="code">norm (x, p) == 1</code>
</p></dd>
</dl>
<p>If <var class="var">A</var> is a vector or a scalar:
</p>
<dl class="table">
<dt><var class="var">p</var> = <code class="code">Inf</code> or <code class="code">"inf"</code></dt>
<dd><p><code class="code">max (abs (<var class="var">A</var>))</code>.
</p>
</dd>
<dt><var class="var">p</var> = <code class="code">-Inf</code></dt>
<dd><p><code class="code">min (abs (<var class="var">A</var>))</code>.
</p>
</dd>
<dt><var class="var">p</var> = <code class="code">"fro"</code></dt>
<dd><p>Frobenius norm of <var class="var">A</var>, <code class="code">sqrt (sumsq (abs (A)))</code>.
</p>
</dd>
<dt><var class="var">p</var> = 0</dt>
<dd><p>Hamming norm—the number of nonzero elements.
</p>
</dd>
<dt>other <var class="var">p</var>, <code class="code"><var class="var">p</var> > 1</code></dt>
<dd><p>p-norm of <var class="var">A</var>, <code class="code">(sum (abs (<var class="var">A</var>) .^ <var class="var">p</var>)) ^ (1/<var class="var">p</var>)</code>.
</p>
</dd>
<dt>other <var class="var">p</var> <code class="code"><var class="var">p</var> < 1</code></dt>
<dd><p>the p-pseudonorm defined as above.
</p></dd>
</dl>
<p>If <var class="var">opt</var> is the value <code class="code">"rows"</code>, treat each row as a vector and
compute its norm. The result is returned as a column vector.
Similarly, if <var class="var">opt</var> is <code class="code">"columns"</code> or <code class="code">"cols"</code> then
compute the norms of each column and return a row vector.
</p>
<p><strong class="strong">See also:</strong> <a class="ref" href="Sparse-Linear-Algebra.html#XREFnormest">normest</a>, <a class="ref" href="Sparse-Linear-Algebra.html#XREFnormest1">normest1</a>, <a class="ref" href="#XREFvecnorm">vecnorm</a>, <a class="ref" href="#XREFcond">cond</a>, <a class="ref" href="Matrix-Factorizations.html#XREFsvd">svd</a>.
</p></dd></dl>
<a class="anchor" id="XREFnull"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-null"><span><code class="def-type"><var class="var">Z</var> =</code> <strong class="def-name">null</strong> <code class="def-code-arguments">(<var class="var">A</var>)</code><a class="copiable-link" href="#index-null"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-null-1"><span><code class="def-type"><var class="var">Z</var> =</code> <strong class="def-name">null</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">tol</var>)</code><a class="copiable-link" href="#index-null-1"> ¶</a></span></dt>
<dd><p>Return an orthonormal basis <var class="var">Z</var> of the null space of <var class="var">A</var>.
</p>
<p>The dimension of the null space <var class="var">Z</var> is taken as the number of singular
values of <var class="var">A</var> not greater than <var class="var">tol</var>. If the argument <var class="var">tol</var>
is missing, it is computed as
</p>
<div class="example">
<pre class="example-preformatted">max (size (<var class="var">A</var>)) * max (svd (<var class="var">A</var>, 0)) * eps
</pre></div>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREForth">orth</a>, <a class="ref" href="Matrix-Factorizations.html#XREFsvd">svd</a>.
</p></dd></dl>
<a class="anchor" id="XREForth"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-orth"><span><code class="def-type"><var class="var">B</var> =</code> <strong class="def-name">orth</strong> <code class="def-code-arguments">(<var class="var">A</var>)</code><a class="copiable-link" href="#index-orth"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-orth-1"><span><code class="def-type"><var class="var">B</var> =</code> <strong class="def-name">orth</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">tol</var>)</code><a class="copiable-link" href="#index-orth-1"> ¶</a></span></dt>
<dd><p>Return an orthonormal basis of the range space of <var class="var">A</var>.
</p>
<p>The dimension of the range space is taken as the number of singular values
of <var class="var">A</var> greater than <var class="var">tol</var>. If the argument <var class="var">tol</var> is missing, it
is computed as
</p>
<div class="example">
<pre class="example-preformatted">max (size (<var class="var">A</var>)) * max (svd (<var class="var">A</var>)) * eps
</pre></div>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFnull">null</a>.
</p></dd></dl>
<a class="anchor" id="XREFmgorth"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-mgorth"><span><code class="def-type">[<var class="var">y</var>, <var class="var">h</var>] =</code> <strong class="def-name">mgorth</strong> <code class="def-code-arguments">(<var class="var">x</var>, <var class="var">v</var>)</code><a class="copiable-link" href="#index-mgorth"> ¶</a></span></dt>
<dd><p>Orthogonalize a given column vector <var class="var">x</var> with respect to a set of
orthonormal vectors comprising the columns of <var class="var">v</var> using the modified
Gram-Schmidt method.
</p>
<p>On exit, <var class="var">y</var> is a unit vector such that:
</p>
<div class="example">
<div class="group"><pre class="example-preformatted"> norm (<var class="var">y</var>) = 1
<var class="var">v</var>' * <var class="var">y</var> = 0
<var class="var">x</var> = [<var class="var">v</var>, <var class="var">y</var>]*<var class="var">h</var>'
</pre></div></div>
</dd></dl>
<a class="anchor" id="XREFpinv"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-pinv"><span><code class="def-type"><var class="var">B</var> =</code> <strong class="def-name">pinv</strong> <code class="def-code-arguments">(<var class="var">A</var>)</code><a class="copiable-link" href="#index-pinv"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-pinv-1"><span><code class="def-type"><var class="var">B</var> =</code> <strong class="def-name">pinv</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">tol</var>)</code><a class="copiable-link" href="#index-pinv-1"> ¶</a></span></dt>
<dd><p>Return the Moore-Penrose pseudoinverse of <var class="var">A</var>.
</p>
<p>Singular values less than <var class="var">tol</var> are ignored.
</p>
<p>If the second argument is omitted, it is taken to be
</p>
<div class="example">
<pre class="example-preformatted">tol = max ([rows(<var class="var">x</var>), columns(<var class="var">x</var>)]) * norm (<var class="var">x</var>) * eps
</pre></div>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFinv">inv</a>, <a class="ref" href="Arithmetic-Ops.html#XREFldivide">ldivide</a>.
</p></dd></dl>
<a class="index-entry-id" id="index-pseudoinverse"></a>
<a class="anchor" id="XREFrank"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-rank"><span><code class="def-type"><var class="var">k</var> =</code> <strong class="def-name">rank</strong> <code class="def-code-arguments">(<var class="var">A</var>)</code><a class="copiable-link" href="#index-rank"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-rank-1"><span><code class="def-type"><var class="var">k</var> =</code> <strong class="def-name">rank</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">tol</var>)</code><a class="copiable-link" href="#index-rank-1"> ¶</a></span></dt>
<dd><p>Compute the rank of matrix <var class="var">A</var>, using the singular value decomposition.
</p>
<p>The rank is taken to be the number of singular values of <var class="var">A</var> that are
greater than the specified tolerance <var class="var">tol</var>. If the second argument is
omitted, it is taken to be
</p>
<div class="example">
<pre class="example-preformatted">tol = max (size (<var class="var">A</var>)) * sigma(1) * eps;
</pre></div>
<p>where <code class="code">eps</code> is machine precision and <code class="code">sigma(1)</code> is the largest
singular value of <var class="var">A</var>.
</p>
<p>The rank of a matrix is the number of linearly independent rows or columns
and equals the dimension of the row and column space. The function
<code class="code">orth</code> may be used to compute an orthonormal basis of the column space.
</p>
<p>For testing if a system <code class="code"><var class="var">A</var>*<var class="var">x</var> = <var class="var">b</var></code> of linear equations
is solvable, one can use
</p>
<div class="example">
<pre class="example-preformatted">rank (<var class="var">A</var>) == rank ([<var class="var">A</var> <var class="var">b</var>])
</pre></div>
<p>In this case, <code class="code"><var class="var">x</var> = <var class="var">A</var> \ <var class="var">b</var></code> finds a particular solution
<var class="var">x</var>. The general solution is <var class="var">x</var> plus the null space of matrix
<var class="var">A</var>. The function <code class="code">null</code> may be used to compute a basis of the
null space.
</p>
<p>Example:
</p>
<div class="example">
<div class="group"><pre class="example-preformatted">A = [1 2 3
4 5 6
7 8 9];
rank (A)
⇒ 2
</pre></div></div>
<p>In this example, the number of linearly independent rows is only 2 because
the final row is a linear combination of the first two rows:
</p>
<div class="example">
<pre class="example-preformatted">A(3,:) == -A(1,:) + 2 * A(2,:)
</pre></div>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFnull">null</a>, <a class="ref" href="#XREForth">orth</a>, <a class="ref" href="Sparse-Linear-Algebra.html#XREFsprank">sprank</a>, <a class="ref" href="Matrix-Factorizations.html#XREFsvd">svd</a>, <a class="ref" href="Mathematical-Constants.html#XREFeps">eps</a>.
</p></dd></dl>
<a class="anchor" id="XREFrcond"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-rcond"><span><code class="def-type"><var class="var">c</var> =</code> <strong class="def-name">rcond</strong> <code class="def-code-arguments">(<var class="var">A</var>)</code><a class="copiable-link" href="#index-rcond"> ¶</a></span></dt>
<dd><p>Compute the 1-norm estimate of the reciprocal condition number as returned
by <small class="sc">LAPACK</small>.
</p>
<p>If the matrix is well-conditioned then <var class="var">c</var> will be near 1 and if the
matrix is poorly conditioned it will be close to 0.
</p>
<p>The matrix <var class="var">A</var> must not be sparse. If the matrix is sparse then
<code class="code">condest (<var class="var">A</var>)</code> or <code class="code">rcond (full (<var class="var">A</var>))</code> should be used
instead.
</p>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFcond">cond</a>, <a class="ref" href="Sparse-Linear-Algebra.html#XREFcondest">condest</a>.
</p></dd></dl>
<a class="anchor" id="XREFtrace"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-trace"><span><code class="def-type"><var class="var">t</var> =</code> <strong class="def-name">trace</strong> <code class="def-code-arguments">(<var class="var">A</var>)</code><a class="copiable-link" href="#index-trace"> ¶</a></span></dt>
<dd><p>Compute the trace of <var class="var">A</var>, the sum of the elements along the main
diagonal.
</p>
<p>The implementation is straightforward: <code class="code">sum (diag (<var class="var">A</var>))</code>.
</p>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFeig">eig</a>.
</p></dd></dl>
<a class="anchor" id="XREFrref"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-rref"><span><code class="def-type"><var class="var">r</var> =</code> <strong class="def-name">rref</strong> <code class="def-code-arguments">(<var class="var">A</var>)</code><a class="copiable-link" href="#index-rref"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-rref-1"><span><code class="def-type"><var class="var">r</var> =</code> <strong class="def-name">rref</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">tol</var>)</code><a class="copiable-link" href="#index-rref-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-rref-2"><span><code class="def-type">[<var class="var">r</var>, <var class="var">k</var>] =</code> <strong class="def-name">rref</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-rref-2"> ¶</a></span></dt>
<dd><p>Return the reduced row echelon form of <var class="var">A</var>.
</p>
<p><var class="var">tol</var> defaults to
<code class="code">eps * max (size (<var class="var">A</var>)) * norm (<var class="var">A</var>, inf)</code>.
</p>
<p>The optional return argument <var class="var">k</var> contains the vector of
"bound variables", which are those columns on which elimination has been
performed.
</p>
</dd></dl>
<a class="anchor" id="XREFvecnorm"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-vecnorm"><span><code class="def-type"><var class="var">n</var> =</code> <strong class="def-name">vecnorm</strong> <code class="def-code-arguments">(<var class="var">A</var>)</code><a class="copiable-link" href="#index-vecnorm"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-vecnorm-1"><span><code class="def-type"><var class="var">n</var> =</code> <strong class="def-name">vecnorm</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">p</var>)</code><a class="copiable-link" href="#index-vecnorm-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-vecnorm-2"><span><code class="def-type"><var class="var">n</var> =</code> <strong class="def-name">vecnorm</strong> <code class="def-code-arguments">(<var class="var">A</var>, <var class="var">p</var>, <var class="var">dim</var>)</code><a class="copiable-link" href="#index-vecnorm-2"> ¶</a></span></dt>
<dd><p>Return the vector p-norm of the elements of array <var class="var">A</var> along dimension
<var class="var">dim</var>.
</p>
<p>The p-norm of a vector is defined as
</p>
<div class="example">
<pre class="example-preformatted"><var class="var">p-norm</var> (<var class="var">A</var>, <var class="var">p</var>) = (sum (abs (<var class="var">A</var>) .^ <var class="var">p</var>)) ^ (1/<var class="var">p</var>)
</pre></div>
<p>The input <var class="var">p</var> must be a positive scalar. If omitted it defaults to 2
(Euclidean norm or distance). Other special values of <var class="var">p</var> are 1
(Manhattan norm, sum of absolute values) and <code class="code">Inf</code> (absolute value of
largest element).
</p>
<p>The input <var class="var">dim</var> specifies the dimension of the array on which the
function operates and must be a positive integer. If omitted the first
non-singleton dimension is used.
</p>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFnorm">norm</a>.
</p></dd></dl>
</div>
<hr>
<div class="nav-panel">
<p>
Next: <a href="Matrix-Factorizations.html">Matrix Factorizations</a>, Previous: <a href="Techniques-Used-for-Linear-Algebra.html">Techniques Used for Linear Algebra</a>, Up: <a href="Linear-Algebra.html">Linear Algebra</a> [<a href="index.html#SEC_Contents" title="Table of contents" rel="contents">Contents</a>][<a href="Concept-Index.html" title="Index" rel="index">Index</a>]</p>
</div>
</body>
</html>
|