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 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209 1210 1211 1212 1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246
|
<HTML>
<HEAD>
<!-- This HTML file has been created by texi2html 1.54+ (gsl)
from ../gsl-ref.texi -->
<TITLE>GNU Scientific Library -- Reference Manual - Linear Algebra</TITLE>
<!-- <LINK rel="stylesheet" title="Default Style Sheet" href="/css/texinfo.css" type="text/css"> -->
<link href="gsl-ref_14.html" rel=Next>
<link href="gsl-ref_12.html" rel=Previous>
<link href="gsl-ref_toc.html" rel=ToC>
</HEAD>
<BODY>
<p>Go to the <A HREF="gsl-ref_1.html">first</A>, <A HREF="gsl-ref_12.html">previous</A>, <A HREF="gsl-ref_14.html">next</A>, <A HREF="gsl-ref_50.html">last</A> section, <A HREF="gsl-ref_toc.html">table of contents</A>.
<P><HR><P>
<H1><A NAME="SEC219" HREF="gsl-ref_toc.html#TOC219">Linear Algebra</A></H1>
<P>
<A NAME="IDX1169"></A>
<A NAME="IDX1170"></A>
<A NAME="IDX1171"></A>
<A NAME="IDX1172"></A>
</P>
<P>
This chapter describes functions for solving linear systems. The
library provides simple linear algebra operations which operate directly
on the <CODE>gsl_vector</CODE> and <CODE>gsl_matrix</CODE> objects. These are
intended for use with "small" systems where simple algorithms are
acceptable.
</P>
<P>
<A NAME="IDX1173"></A>
Anyone interested in large systems will want to use the sophisticated
routines found in LAPACK. The Fortran version of LAPACK is
recommended as the standard package for linear algebra. It supports
blocked algorithms, specialized data representations and other
optimizations.
</P>
<P>
The functions described in this chapter are declared in the header file
<TT>'gsl_linalg.h'</TT>.
</P>
<H2><A NAME="SEC220" HREF="gsl-ref_toc.html#TOC220">LU Decomposition</A></H2>
<P>
<A NAME="IDX1174"></A>
</P>
<P>
A general square matrix A has an LU decomposition into
upper and lower triangular matrices,
</P>
<PRE class="example">
P A = L U
</PRE>
<P>
where P is a permutation matrix, L is unit lower
triangular matrix and U is upper triangular matrix. For square
matrices this decomposition can be used to convert the linear system
A x = b into a pair of triangular systems (L y = P b,
U x = y), which can be solved by forward and back-substitution.
Note that the LU decomposition is valid for singular matrices.
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_LU_decomp</B> <I>(gsl_matrix * <VAR>A</VAR>, gsl_permutation * <VAR>p</VAR>, int *<VAR>signum</VAR>)</I>
<DD><A NAME="IDX1175"></A>
<DT><U>Function:</U> int <B>gsl_linalg_complex_LU_decomp</B> <I>(gsl_matrix_complex * <VAR>A</VAR>, gsl_permutation * <VAR>p</VAR>, int *<VAR>signum</VAR>)</I>
<DD><A NAME="IDX1176"></A>
These functions factorize the square matrix <VAR>A</VAR> into the LU
decomposition PA = LU. On output the diagonal and upper
triangular part of the input matrix <VAR>A</VAR> contain the matrix
U. The lower triangular part of the input matrix (excluding the
diagonal) contains L. The diagonal elements of L are
unity, and are not stored.
</P>
<P>
The permutation matrix P is encoded in the permutation
<VAR>p</VAR>. The j-th column of the matrix P is given by the
k-th column of the identity matrix, where k = p_j the
j-th element of the permutation vector. The sign of the
permutation is given by <VAR>signum</VAR>. It has the value (-1)^n,
where n is the number of interchanges in the permutation.
</P>
<P>
The algorithm used in the decomposition is Gaussian Elimination with
partial pivoting (Golub & Van Loan, <CITE>Matrix Computations</CITE>,
Algorithm 3.4.1).
</DL>
</P>
<P>
<A NAME="IDX1177"></A>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_LU_solve</B> <I>(const gsl_matrix * <VAR>LU</VAR>, const gsl_permutation * <VAR>p</VAR>, const gsl_vector * <VAR>b</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1178"></A>
<DT><U>Function:</U> int <B>gsl_linalg_complex_LU_solve</B> <I>(const gsl_matrix_complex * <VAR>LU</VAR>, const gsl_permutation * <VAR>p</VAR>, const gsl_vector_complex * <VAR>b</VAR>, gsl_vector_complex * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1179"></A>
These functions solve the system A x = b using the LU
decomposition of A into (<VAR>LU</VAR>, <VAR>p</VAR>) given by
<CODE>gsl_linalg_LU_decomp</CODE> or <CODE>gsl_linalg_complex_LU_decomp</CODE>.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_LU_svx</B> <I>(const gsl_matrix * <VAR>LU</VAR>, const gsl_permutation * <VAR>p</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1180"></A>
<DT><U>Function:</U> int <B>gsl_linalg_complex_LU_svx</B> <I>(const gsl_matrix_complex * <VAR>LU</VAR>, const gsl_permutation * <VAR>p</VAR>, gsl_vector_complex * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1181"></A>
These functions solve the system A x = b in-place using the
LU decomposition of A into (<VAR>LU</VAR>,<VAR>p</VAR>). On input
<VAR>x</VAR> should contain the right-hand side b, which is replaced
by the solution on output.
</DL>
</P>
<P>
<A NAME="IDX1182"></A>
<A NAME="IDX1183"></A>
<A NAME="IDX1184"></A>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_LU_refine</B> <I>(const gsl_matrix * <VAR>A</VAR>, const gsl_matrix * <VAR>LU</VAR>, const gsl_permutation * <VAR>p</VAR>, const gsl_vector * <VAR>b</VAR>, gsl_vector * <VAR>x</VAR>, gsl_vector * <VAR>residual</VAR>)</I>
<DD><A NAME="IDX1185"></A>
<DT><U>Function:</U> int <B>gsl_linalg_complex_LU_refine</B> <I>(const gsl_matrix_complex * <VAR>A</VAR>, const gsl_matrix_complex * <VAR>LU</VAR>, const gsl_permutation * <VAR>p</VAR>, const gsl_vector_complex * <VAR>b</VAR>, gsl_vector_complex * <VAR>x</VAR>, gsl_vector_complex * <VAR>residual</VAR>)</I>
<DD><A NAME="IDX1186"></A>
These functions apply an iterative improvement to <VAR>x</VAR>, the solution
of A x = b, using the LU decomposition of A into
(<VAR>LU</VAR>,<VAR>p</VAR>). The initial residual r = A x - b is also
computed and stored in <VAR>residual</VAR>.
</DL>
</P>
<P>
<A NAME="IDX1187"></A>
<A NAME="IDX1188"></A>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_LU_invert</B> <I>(const gsl_matrix * <VAR>LU</VAR>, const gsl_permutation * <VAR>p</VAR>, gsl_matrix * <VAR>inverse</VAR>)</I>
<DD><A NAME="IDX1189"></A>
<DT><U>Function:</U> int <B>gsl_linalg_complex_LU_invert</B> <I>(const gsl_matrix_complex * <VAR>LU</VAR>, const gsl_permutation * <VAR>p</VAR>, gsl_matrix_complex * <VAR>inverse</VAR>)</I>
<DD><A NAME="IDX1190"></A>
These functions compute the inverse of a matrix A from its
LU decomposition (<VAR>LU</VAR>,<VAR>p</VAR>), storing the result in the
matrix <VAR>inverse</VAR>. The inverse is computed by solving the system
A x = b for each column of the identity matrix. It is preferable
to avoid direct computation of the inverse whenever possible.
</DL>
</P>
<P>
<A NAME="IDX1191"></A>
<A NAME="IDX1192"></A>
<DL>
<DT><U>Function:</U> double <B>gsl_linalg_LU_det</B> <I>(gsl_matrix * <VAR>LU</VAR>, int <VAR>signum</VAR>)</I>
<DD><A NAME="IDX1193"></A>
<DT><U>Function:</U> gsl_complex <B>gsl_linalg_complex_LU_det</B> <I>(gsl_matrix_complex * <VAR>LU</VAR>, int <VAR>signum</VAR>)</I>
<DD><A NAME="IDX1194"></A>
These functions compute the determinant of a matrix A from its
LU decomposition, <VAR>LU</VAR>. The determinant is computed as the
product of the diagonal elements of U and the sign of the row
permutation <VAR>signum</VAR>.
</DL>
</P>
<P>
<A NAME="IDX1195"></A>
<DL>
<DT><U>Function:</U> double <B>gsl_linalg_LU_lndet</B> <I>(gsl_matrix * <VAR>LU</VAR>)</I>
<DD><A NAME="IDX1196"></A>
<DT><U>Function:</U> double <B>gsl_linalg_complex_LU_lndet</B> <I>(gsl_matrix_complex * <VAR>LU</VAR>)</I>
<DD><A NAME="IDX1197"></A>
These functions compute the logarithm of the absolute value of the
determinant of a matrix A, \ln|det(A)|, from its LU
decomposition, <VAR>LU</VAR>. This function may be useful if the direct
computation of the determinant would overflow or underflow.
</DL>
</P>
<P>
<A NAME="IDX1198"></A>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_LU_sgndet</B> <I>(gsl_matrix * <VAR>LU</VAR>, int <VAR>signum</VAR>)</I>
<DD><A NAME="IDX1199"></A>
<DT><U>Function:</U> gsl_complex <B>gsl_linalg_complex_LU_sgndet</B> <I>(gsl_matrix_complex * <VAR>LU</VAR>, int <VAR>signum</VAR>)</I>
<DD><A NAME="IDX1200"></A>
These functions compute the sign or phase factor of the determinant of a
matrix A, det(A)/|det(A)|, from its LU decomposition,
<VAR>LU</VAR>.
</DL>
</P>
<H2><A NAME="SEC221" HREF="gsl-ref_toc.html#TOC221">QR Decomposition</A></H2>
<P>
<A NAME="IDX1201"></A>
</P>
<P>
A general rectangular M-by-N matrix A has a
QR decomposition into the product of an orthogonal
M-by-M square matrix Q (where Q^T Q = I) and
an M-by-N right-triangular matrix R,
</P>
<PRE class="example">
A = Q R
</PRE>
<P>
This decomposition can be used to convert the linear system A x =
b into the triangular system R x = Q^T b, which can be solved by
back-substitution. Another use of the QR decomposition is to
compute an orthonormal basis for a set of vectors. The first N
columns of Q form an orthonormal basis for the range of A,
ran(A), when A has full column rank.
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_QR_decomp</B> <I>(gsl_matrix * <VAR>A</VAR>, gsl_vector * <VAR>tau</VAR>)</I>
<DD><A NAME="IDX1202"></A>
This function factorizes the M-by-N matrix <VAR>A</VAR> into
the QR decomposition A = Q R. On output the diagonal and
upper triangular part of the input matrix contain the matrix
R. The vector <VAR>tau</VAR> and the columns of the lower triangular
part of the matrix <VAR>A</VAR> contain the Householder coefficients and
Householder vectors which encode the orthogonal matrix <VAR>Q</VAR>. The
vector <VAR>tau</VAR> must be of length k=\min(M,N). The matrix
Q is related to these components by, Q = Q_k ... Q_2 Q_1
where Q_i = I - \tau_i v_i v_i^T and v_i is the
Householder vector v_i =
(0,...,1,A(i+1,i),A(i+2,i),...,A(m,i)). This is the same storage scheme
as used by LAPACK.
</P>
<P>
The algorithm used to perform the decomposition is Householder QR (Golub
& Van Loan, <CITE>Matrix Computations</CITE>, Algorithm 5.2.1).
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_QR_solve</B> <I>(const gsl_matrix * <VAR>QR</VAR>, const gsl_vector * <VAR>tau</VAR>, const gsl_vector * <VAR>b</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1203"></A>
This function solves the system A x = b using the QR
decomposition of A into (<VAR>QR</VAR>, <VAR>tau</VAR>) given by
<CODE>gsl_linalg_QR_decomp</CODE>.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_QR_svx</B> <I>(const gsl_matrix * <VAR>QR</VAR>, const gsl_vector * <VAR>tau</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1204"></A>
This function solves the system A x = b in-place using the
QR decomposition of A into (<VAR>QR</VAR>,<VAR>tau</VAR>) given by
<CODE>gsl_linalg_QR_decomp</CODE>. On input <VAR>x</VAR> should contain the
right-hand side b, which is replaced by the solution on output.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_QR_lssolve</B> <I>(const gsl_matrix * <VAR>QR</VAR>, const gsl_vector * <VAR>tau</VAR>, const gsl_vector * <VAR>b</VAR>, gsl_vector * <VAR>x</VAR>, gsl_vector * <VAR>residual</VAR>)</I>
<DD><A NAME="IDX1205"></A>
This function finds the least squares solution to the overdetermined
system A x = b where the matrix <VAR>A</VAR> has more rows than
columns. The least squares solution minimizes the Euclidean norm of the
residual, ||Ax - b||.The routine uses the QR decomposition
of A into (<VAR>QR</VAR>, <VAR>tau</VAR>) given by
<CODE>gsl_linalg_QR_decomp</CODE>. The solution is returned in <VAR>x</VAR>. The
residual is computed as a by-product and stored in <VAR>residual</VAR>.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_QR_QTvec</B> <I>(const gsl_matrix * <VAR>QR</VAR>, const gsl_vector * <VAR>tau</VAR>, gsl_vector * <VAR>v</VAR>)</I>
<DD><A NAME="IDX1206"></A>
This function applies the matrix Q^T encoded in the decomposition
(<VAR>QR</VAR>,<VAR>tau</VAR>) to the vector <VAR>v</VAR>, storing the result Q^T
v in <VAR>v</VAR>. The matrix multiplication is carried out directly using
the encoding of the Householder vectors without needing to form the full
matrix Q^T.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_QR_Qvec</B> <I>(const gsl_matrix * <VAR>QR</VAR>, const gsl_vector * <VAR>tau</VAR>, gsl_vector * <VAR>v</VAR>)</I>
<DD><A NAME="IDX1207"></A>
This function applies the matrix Q encoded in the decomposition
(<VAR>QR</VAR>,<VAR>tau</VAR>) to the vector <VAR>v</VAR>, storing the result Q
v in <VAR>v</VAR>. The matrix multiplication is carried out directly using
the encoding of the Householder vectors without needing to form the full
matrix Q.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_QR_Rsolve</B> <I>(const gsl_matrix * <VAR>QR</VAR>, const gsl_vector * <VAR>b</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1208"></A>
This function solves the triangular system R x = b for
<VAR>x</VAR>. It may be useful if the product b' = Q^T b has already
been computed using <CODE>gsl_linalg_QR_QTvec</CODE>.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_QR_Rsvx</B> <I>(const gsl_matrix * <VAR>QR</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1209"></A>
This function solves the triangular system R x = b for <VAR>x</VAR>
in-place. On input <VAR>x</VAR> should contain the right-hand side b
and is replaced by the solution on output. This function may be useful if
the product b' = Q^T b has already been computed using
<CODE>gsl_linalg_QR_QTvec</CODE>.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_QR_unpack</B> <I>(const gsl_matrix * <VAR>QR</VAR>, const gsl_vector * <VAR>tau</VAR>, gsl_matrix * <VAR>Q</VAR>, gsl_matrix * <VAR>R</VAR>)</I>
<DD><A NAME="IDX1210"></A>
This function unpacks the encoded QR decomposition
(<VAR>QR</VAR>,<VAR>tau</VAR>) into the matrices <VAR>Q</VAR> and <VAR>R</VAR>, where
<VAR>Q</VAR> is M-by-M and <VAR>R</VAR> is M-by-N.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_QR_QRsolve</B> <I>(gsl_matrix * <VAR>Q</VAR>, gsl_matrix * <VAR>R</VAR>, const gsl_vector * <VAR>b</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1211"></A>
This function solves the system R x = Q^T b for <VAR>x</VAR>. It can
be used when the QR decomposition of a matrix is available in
unpacked form as (<VAR>Q</VAR>, <VAR>R</VAR>).
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_QR_update</B> <I>(gsl_matrix * <VAR>Q</VAR>, gsl_matrix * <VAR>R</VAR>, gsl_vector * <VAR>w</VAR>, const gsl_vector * <VAR>v</VAR>)</I>
<DD><A NAME="IDX1212"></A>
This function performs a rank-1 update w v^T of the QR
decomposition (<VAR>Q</VAR>, <VAR>R</VAR>). The update is given by Q'R' = Q
R + w v^T where the output matrices Q' and R' are also
orthogonal and right triangular. Note that <VAR>w</VAR> is destroyed by the
update.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_R_solve</B> <I>(const gsl_matrix * <VAR>R</VAR>, const gsl_vector * <VAR>b</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1213"></A>
This function solves the triangular system R x = b for the
N-by-N matrix <VAR>R</VAR>.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_R_svx</B> <I>(const gsl_matrix * <VAR>R</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1214"></A>
This function solves the triangular system R x = b in-place. On
input <VAR>x</VAR> should contain the right-hand side b, which is
replaced by the solution on output.
</DL>
</P>
<H2><A NAME="SEC222" HREF="gsl-ref_toc.html#TOC222">QR Decomposition with Column Pivoting</A></H2>
<P>
<A NAME="IDX1215"></A>
</P>
<P>
The QR decomposition can be extended to the rank deficient case
by introducing a column permutation P,
</P>
<PRE class="example">
A P = Q R
</PRE>
<P>
The first r columns of this Q form an orthonormal basis
for the range of A for a matrix with column rank r. This
decomposition can also be used to convert the linear system A x =
b into the triangular system R y = Q^T b, x = P y, which can be
solved by back-substitution and permutation. We denote the QR
decomposition with column pivoting by QRP^T since A = Q R
P^T.
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_QRPT_decomp</B> <I>(gsl_matrix * <VAR>A</VAR>, gsl_vector * <VAR>tau</VAR>, gsl_permutation * <VAR>p</VAR>, int *<VAR>signum</VAR>, gsl_vector * <VAR>norm</VAR>)</I>
<DD><A NAME="IDX1216"></A>
This function factorizes the M-by-N matrix <VAR>A</VAR> into
the QRP^T decomposition A = Q R P^T. On output the
diagonal and upper triangular part of the input matrix contain the
matrix R. The permutation matrix P is stored in the
permutation <VAR>p</VAR>. The sign of the permutation is given by
<VAR>signum</VAR>. It has the value (-1)^n, where n is the
number of interchanges in the permutation. The vector <VAR>tau</VAR> and the
columns of the lower triangular part of the matrix <VAR>A</VAR> contain the
Householder coefficients and vectors which encode the orthogonal matrix
<VAR>Q</VAR>. The vector <VAR>tau</VAR> must be of length k=\min(M,N). The
matrix Q is related to these components by, Q = Q_k ... Q_2
Q_1 where Q_i = I - \tau_i v_i v_i^T and v_i is the
Householder vector v_i =
(0,...,1,A(i+1,i),A(i+2,i),...,A(m,i)). This is the same storage scheme
as used by LAPACK. The vector <VAR>norm</VAR> is workspace of length
<VAR>N</VAR> used for column pivoting.
</P>
<P>
The algorithm used to perform the decomposition is Householder QR with
column pivoting (Golub & Van Loan, <CITE>Matrix Computations</CITE>, Algorithm
5.4.1).
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_QRPT_decomp2</B> <I>(const gsl_matrix * <VAR>A</VAR>, gsl_matrix * <VAR>q</VAR>, gsl_matrix * <VAR>r</VAR>, gsl_vector * <VAR>tau</VAR>, gsl_permutation * <VAR>p</VAR>, int *<VAR>signum</VAR>, gsl_vector * <VAR>norm</VAR>)</I>
<DD><A NAME="IDX1217"></A>
This function factorizes the matrix <VAR>A</VAR> into the decomposition
A = Q R P^T without modifying <VAR>A</VAR> itself and storing the
output in the separate matrices <VAR>q</VAR> and <VAR>r</VAR>.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_QRPT_solve</B> <I>(const gsl_matrix * <VAR>QR</VAR>, const gsl_vector * <VAR>tau</VAR>, const gsl_permutation * <VAR>p</VAR>, const gsl_vector * <VAR>b</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1218"></A>
This function solves the system A x = b using the QRP^T
decomposition of A into (<VAR>QR</VAR>, <VAR>tau</VAR>, <VAR>p</VAR>) given by
<CODE>gsl_linalg_QRPT_decomp</CODE>.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_QRPT_svx</B> <I>(const gsl_matrix * <VAR>QR</VAR>, const gsl_vector * <VAR>tau</VAR>, const gsl_permutation * <VAR>p</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1219"></A>
This function solves the system A x = b in-place using the
QRP^T decomposition of A into
(<VAR>QR</VAR>,<VAR>tau</VAR>,<VAR>p</VAR>). On input <VAR>x</VAR> should contain the
right-hand side b, which is replaced by the solution on output.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_QRPT_QRsolve</B> <I>(const gsl_matrix * <VAR>Q</VAR>, const gsl_matrix * <VAR>R</VAR>, const gsl_permutation * <VAR>p</VAR>, const gsl_vector * <VAR>b</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1220"></A>
This function solves the system R P^T x = Q^T b for <VAR>x</VAR>. It can
be used when the QR decomposition of a matrix is available in
unpacked form as (<VAR>Q</VAR>, <VAR>R</VAR>).
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_QRPT_update</B> <I>(gsl_matrix * <VAR>Q</VAR>, gsl_matrix * <VAR>R</VAR>, const gsl_permutation * <VAR>p</VAR>, gsl_vector * <VAR>u</VAR>, const gsl_vector * <VAR>v</VAR>)</I>
<DD><A NAME="IDX1221"></A>
This function performs a rank-1 update w v^T of the QRP^T
decomposition (<VAR>Q</VAR>, <VAR>R</VAR>, <VAR>p</VAR>). The update is given by
Q'R' = Q R + w v^T where the output matrices Q' and
R' are also orthogonal and right triangular. Note that <VAR>w</VAR> is
destroyed by the update. The permutation <VAR>p</VAR> is not changed.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_QRPT_Rsolve</B> <I>(const gsl_matrix * <VAR>QR</VAR>, const gsl_permutation * <VAR>p</VAR>, const gsl_vector * <VAR>b</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1222"></A>
This function solves the triangular system R P^T x = b for the
N-by-N matrix R contained in <VAR>QR</VAR>.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_QRPT_Rsvx</B> <I>(const gsl_matrix * <VAR>QR</VAR>, const gsl_permutation * <VAR>p</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1223"></A>
This function solves the triangular system R P^T x = b in-place
for the N-by-N matrix R contained in <VAR>QR</VAR>. On
input <VAR>x</VAR> should contain the right-hand side b, which is
replaced by the solution on output.
</DL>
</P>
<H2><A NAME="SEC223" HREF="gsl-ref_toc.html#TOC223">Singular Value Decomposition</A></H2>
<P>
<A NAME="IDX1224"></A>
<A NAME="IDX1225"></A>
</P>
<P>
A general rectangular M-by-N matrix A has a
singular value decomposition (SVD) into the product of an
M-by-N orthogonal matrix U, an N-by-N
diagonal matrix of singular values S and the transpose of an
N-by-N orthogonal square matrix V,
</P>
<PRE class="example">
A = U S V^T
</PRE>
<P>
The singular values
\sigma_i = S_{ii} are all non-negative and are
generally chosen to form a non-increasing sequence
\sigma_1 >= \sigma_2 >= ... >= \sigma_N >= 0.
</P>
<P>
The singular value decomposition of a matrix has many practical uses.
The condition number of the matrix is given by the ratio of the largest
singular value to the smallest singular value. The presence of a zero
singular value indicates that the matrix is singular. The number of
non-zero singular values indicates the rank of the matrix. In practice
singular value decomposition of a rank-deficient matrix will not produce
exact zeroes for singular values, due to finite numerical
precision. Small singular values should be edited by choosing a suitable
tolerance.
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_SV_decomp</B> <I>(gsl_matrix * <VAR>A</VAR>, gsl_matrix * <VAR>V</VAR>, gsl_vector * <VAR>S</VAR>, gsl_vector * <VAR>work</VAR>)</I>
<DD><A NAME="IDX1226"></A>
This function factorizes the M-by-N matrix <VAR>A</VAR> into
the singular value decomposition A = U S V^T. On output the
matrix <VAR>A</VAR> is replaced by U. The diagonal elements of the
singular value matrix S are stored in the vector <VAR>S</VAR>. The
singular values are non-negative and form a non-increasing sequence from
S_1 to S_N. The matrix <VAR>V</VAR> contains the elements of
V in untransposed form. To form the product U S V^T it is
necessary to take the transpose of <VAR>V</VAR>. A workspace of length
<VAR>N</VAR> is required in <VAR>work</VAR>.
</P>
<P>
This routine uses the Golub-Reinsch SVD algorithm.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_SV_decomp_mod</B> <I>(gsl_matrix * <VAR>A</VAR>, gsl_matrix * <VAR>X</VAR>, gsl_matrix * <VAR>V</VAR>, gsl_vector * <VAR>S</VAR>, gsl_vector * <VAR>work</VAR>)</I>
<DD><A NAME="IDX1227"></A>
This function computes the SVD using the modified Golub-Reinsch
algorithm, which is faster for
M>>N. It requires the vector
<VAR>work</VAR> and the N-by-N matrix <VAR>X</VAR> as additional
working space.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_SV_decomp_jacobi</B> <I>(gsl_matrix * <VAR>A</VAR>, gsl_matrix * <VAR>V</VAR>, gsl_vector * <VAR>S</VAR>)</I>
<DD><A NAME="IDX1228"></A>
This function computes the SVD using one-sided Jacobi orthogonalization
(see references for details). The Jacobi method can compute singular
values to higher relative accuracy than Golub-Reinsch algorithms.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_SV_solve</B> <I>(gsl_matrix * <VAR>U</VAR>, gsl_matrix * <VAR>V</VAR>, gsl_vector * <VAR>S</VAR>, const gsl_vector * <VAR>b</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1229"></A>
This function solves the system A x = b using the singular value
decomposition (<VAR>U</VAR>, <VAR>S</VAR>, <VAR>V</VAR>) of A given by
<CODE>gsl_linalg_SV_decomp</CODE>.
</P>
<P>
Only non-zero singular values are used in computing the solution. The
parts of the solution corresponding to singular values of zero are
ignored. Other singular values can be edited out by setting them to
zero before calling this function.
</P>
<P>
In the over-determined case where <VAR>A</VAR> has more rows than columns the
system is solved in the least squares sense, returning the solution
<VAR>x</VAR> which minimizes ||A x - b||_2.
</DL>
</P>
<H2><A NAME="SEC224" HREF="gsl-ref_toc.html#TOC224">Cholesky Decomposition</A></H2>
<P>
<A NAME="IDX1230"></A>
<A NAME="IDX1231"></A>
<A NAME="IDX1232"></A>
</P>
<P>
A symmetric, positive definite square matrix A has a Cholesky
decomposition into a product of a lower triangular matrix L and
its transpose L^T,
</P>
<PRE class="example">
A = L L^T
</PRE>
<P>
This is sometimes referred to as taking the square-root of a matrix. The
Cholesky decomposition can only be carried out when all the eigenvalues
of the matrix are positive. This decomposition can be used to convert
the linear system A x = b into a pair of triangular systems
(L y = b, L^T x = y), which can be solved by forward and
back-substitution.
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_cholesky_decomp</B> <I>(gsl_matrix * <VAR>A</VAR>)</I>
<DD><A NAME="IDX1233"></A>
This function factorizes the positive-definite square matrix <VAR>A</VAR>
into the Cholesky decomposition A = L L^T. On output the diagonal
and lower triangular part of the input matrix <VAR>A</VAR> contain the matrix
L. The upper triangular part of the input matrix contains
L^T, the diagonal terms being identical for both L and
L^T. If the matrix is not positive-definite then the
decomposition will fail, returning the error code <CODE>GSL_EDOM</CODE>.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_cholesky_solve</B> <I>(const gsl_matrix * <VAR>cholesky</VAR>, const gsl_vector * <VAR>b</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1234"></A>
This function solves the system A x = b using the Cholesky
decomposition of A into the matrix <VAR>cholesky</VAR> given by
<CODE>gsl_linalg_cholesky_decomp</CODE>.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_cholesky_svx</B> <I>(const gsl_matrix * <VAR>cholesky</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1235"></A>
This function solves the system A x = b in-place using the
Cholesky decomposition of A into the matrix <VAR>cholesky</VAR> given
by <CODE>gsl_linalg_cholesky_decomp</CODE>. On input <VAR>x</VAR> should contain
the right-hand side b, which is replaced by the solution on
output.
</DL>
</P>
<H2><A NAME="SEC225" HREF="gsl-ref_toc.html#TOC225">Tridiagonal Decomposition of Real Symmetric Matrices</A></H2>
<P>
<A NAME="IDX1236"></A>
</P>
<P>
A symmetric matrix A can be factorized by similarity
transformations into the form,
</P>
<PRE class="example">
A = Q T Q^T
</PRE>
<P>
where Q is an orthogonal matrix and T is a symmetric
tridiagonal matrix.
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_symmtd_decomp</B> <I>(gsl_matrix * <VAR>A</VAR>, gsl_vector * <VAR>tau</VAR>)</I>
<DD><A NAME="IDX1237"></A>
This function factorizes the symmetric square matrix <VAR>A</VAR> into the
symmetric tridiagonal decomposition Q T Q^T. On output the
diagonal and subdiagonal part of the input matrix <VAR>A</VAR> contain the
tridiagonal matrix T. The remaining lower triangular part of the
input matrix contains the Householder vectors which, together with the
Householder coefficients <VAR>tau</VAR>, encode the orthogonal matrix
Q. This storage scheme is the same as used by LAPACK. The
upper triangular part of <VAR>A</VAR> is not referenced.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_symmtd_unpack</B> <I>(const gsl_matrix * <VAR>A</VAR>, const gsl_vector * <VAR>tau</VAR>, gsl_matrix * <VAR>Q</VAR>, gsl_vector * <VAR>diag</VAR>, gsl_vector * <VAR>subdiag</VAR>)</I>
<DD><A NAME="IDX1238"></A>
This function unpacks the encoded symmetric tridiagonal decomposition
(<VAR>A</VAR>, <VAR>tau</VAR>) obtained from <CODE>gsl_linalg_symmtd_decomp</CODE> into
the orthogonal matrix <VAR>Q</VAR>, the vector of diagonal elements <VAR>diag</VAR>
and the vector of subdiagonal elements <VAR>subdiag</VAR>.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_symmtd_unpack_T</B> <I>(const gsl_matrix * <VAR>A</VAR>, gsl_vector * <VAR>diag</VAR>, gsl_vector * <VAR>subdiag</VAR>)</I>
<DD><A NAME="IDX1239"></A>
This function unpacks the diagonal and subdiagonal of the encoded
symmetric tridiagonal decomposition (<VAR>A</VAR>, <VAR>tau</VAR>) obtained from
<CODE>gsl_linalg_symmtd_decomp</CODE> into the vectors <VAR>diag</VAR> and <VAR>subdiag</VAR>.
</DL>
</P>
<H2><A NAME="SEC226" HREF="gsl-ref_toc.html#TOC226">Tridiagonal Decomposition of Hermitian Matrices</A></H2>
<P>
<A NAME="IDX1240"></A>
</P>
<P>
A hermitian matrix A can be factorized by similarity
transformations into the form,
</P>
<PRE class="example">
A = U T U^T
</PRE>
<P>
where U is an unitary matrix and T is a real symmetric
tridiagonal matrix.
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_hermtd_decomp</B> <I>(gsl_matrix_complex * <VAR>A</VAR>, gsl_vector_complex * <VAR>tau</VAR>)</I>
<DD><A NAME="IDX1241"></A>
This function factorizes the hermitian matrix <VAR>A</VAR> into the symmetric
tridiagonal decomposition U T U^T. On output the real parts of
the diagonal and subdiagonal part of the input matrix <VAR>A</VAR> contain
the tridiagonal matrix T. The remaining lower triangular part of
the input matrix contains the Householder vectors which, together with
the Householder coefficients <VAR>tau</VAR>, encode the orthogonal matrix
Q. This storage scheme is the same as used by LAPACK. The
upper triangular part of <VAR>A</VAR> and imaginary parts of the diagonal are
not referenced.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_hermtd_unpack</B> <I>(const gsl_matrix_complex * <VAR>A</VAR>, const gsl_vector_complex * <VAR>tau</VAR>, gsl_matrix_complex * <VAR>Q</VAR>, gsl_vector * <VAR>diag</VAR>, gsl_vector * <VAR>subdiag</VAR>)</I>
<DD><A NAME="IDX1242"></A>
This function unpacks the encoded tridiagonal decomposition (<VAR>A</VAR>,
<VAR>tau</VAR>) obtained from <CODE>gsl_linalg_hermtd_decomp</CODE> into the
unitary matrix <VAR>U</VAR>, the real vector of diagonal elements <VAR>diag</VAR> and
the real vector of subdiagonal elements <VAR>subdiag</VAR>.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_hermtd_unpack_T</B> <I>(const gsl_matrix_complex * <VAR>A</VAR>, gsl_vector * <VAR>diag</VAR>, gsl_vector * <VAR>subdiag</VAR>)</I>
<DD><A NAME="IDX1243"></A>
This function unpacks the diagonal and subdiagonal of the encoded
tridiagonal decomposition (<VAR>A</VAR>, <VAR>tau</VAR>) obtained from
<CODE>gsl_linalg_hermtd_decomp</CODE> into the real vectors <VAR>diag</VAR> and
<VAR>subdiag</VAR>.
</DL>
</P>
<H2><A NAME="SEC227" HREF="gsl-ref_toc.html#TOC227">Bidiagonalization</A></H2>
<P>
<A NAME="IDX1244"></A>
</P>
<P>
A general matrix A can be factorized by similarity
transformations into the form,
</P>
<PRE class="example">
A = U B V^T
</PRE>
<P>
where U and V are orthogonal matrices and B is a
N-by-N bidiagonal matrix with non-zero entries only on the
diagonal and superdiagonal. The size of <VAR>U</VAR> is M-by-N
and the size of <VAR>V</VAR> is N-by-N.
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_bidiag_decomp</B> <I>(gsl_matrix * <VAR>A</VAR>, gsl_vector * <VAR>tau_U</VAR>, gsl_vector * <VAR>tau_V</VAR>)</I>
<DD><A NAME="IDX1245"></A>
This function factorizes the M-by-N matrix <VAR>A</VAR> into
bidiagonal form U B V^T. The diagonal and superdiagonal of the
matrix B are stored in the diagonal and superdiagonal of <VAR>A</VAR>.
The orthogonal matrices U and <VAR>V</VAR> are stored as compressed
Householder vectors in the remaining elements of <VAR>A</VAR>. The
Householder coefficients are stored in the vectors <VAR>tau_U</VAR> and
<VAR>tau_V</VAR>. The length of <VAR>tau_U</VAR> must equal the number of
elements in the diagonal of <VAR>A</VAR> and the length of <VAR>tau_V</VAR> should
be one element shorter.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_bidiag_unpack</B> <I>(const gsl_matrix * <VAR>A</VAR>, const gsl_vector * <VAR>tau_U</VAR>, gsl_matrix * <VAR>U</VAR>, const gsl_vector * <VAR>tau_V</VAR>, gsl_matrix * <VAR>V</VAR>, gsl_vector * <VAR>diag</VAR>, gsl_vector * <VAR>superdiag</VAR>)</I>
<DD><A NAME="IDX1246"></A>
This function unpacks the bidiagonal decomposition of <VAR>A</VAR> given by
<CODE>gsl_linalg_bidiag_decomp</CODE>, (<VAR>A</VAR>, <VAR>tau_U</VAR>, <VAR>tau_V</VAR>)
into the separate orthogonal matrices <VAR>U</VAR>, <VAR>V</VAR> and the diagonal
vector <VAR>diag</VAR> and superdiagonal <VAR>superdiag</VAR>. Note that <VAR>U</VAR>
is stored as a compact M-by-N orthogonal matrix satisfying
U^T U = I for efficiency.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_bidiag_unpack2</B> <I>(gsl_matrix * <VAR>A</VAR>, gsl_vector * <VAR>tau_U</VAR>, gsl_vector * <VAR>tau_V</VAR>, gsl_matrix * <VAR>V</VAR>)</I>
<DD><A NAME="IDX1247"></A>
This function unpacks the bidiagonal decomposition of <VAR>A</VAR> given by
<CODE>gsl_linalg_bidiag_decomp</CODE>, (<VAR>A</VAR>, <VAR>tau_U</VAR>, <VAR>tau_V</VAR>)
into the separate orthogonal matrices <VAR>U</VAR>, <VAR>V</VAR> and the diagonal
vector <VAR>diag</VAR> and superdiagonal <VAR>superdiag</VAR>. The matrix <VAR>U</VAR>
is stored in-place in <VAR>A</VAR>.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_bidiag_unpack_B</B> <I>(const gsl_matrix * <VAR>A</VAR>, gsl_vector * <VAR>diag</VAR>, gsl_vector * <VAR>superdiag</VAR>)</I>
<DD><A NAME="IDX1248"></A>
This function unpacks the diagonal and superdiagonal of the bidiagonal
decomposition of <VAR>A</VAR> given by <CODE>gsl_linalg_bidiag_decomp</CODE>, into
the diagonal vector <VAR>diag</VAR> and superdiagonal vector <VAR>superdiag</VAR>.
</DL>
</P>
<H2><A NAME="SEC228" HREF="gsl-ref_toc.html#TOC228">Householder Transformations</A></H2>
<P>
<A NAME="IDX1249"></A>
<A NAME="IDX1250"></A>
<A NAME="IDX1251"></A>
</P>
<P>
A Householder transformation is a rank-1 modification of the identity
matrix which can be used to zero out selected elements of a vector. A
Householder matrix P takes the form,
</P>
<PRE class="example">
P = I - \tau v v^T
</PRE>
<P>
where v is a vector (called the <I>Householder vector</I>) and
\tau = 2/(v^T v). The functions described in this section use the
rank-1 structure of the Householder matrix to create and apply
Householder transformations efficiently.
</P>
<P>
<DL>
<DT><U>Function:</U> double <B>gsl_linalg_householder_transform</B> <I>(gsl_vector * <VAR>v</VAR>)</I>
<DD><A NAME="IDX1252"></A>
This function prepares a Householder transformation P = I - \tau v
v^T which can be used to zero all the elements of the input vector except
the first. On output the transformation is stored in the vector <VAR>v</VAR>
and the scalar \tau is returned.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_householder_hm</B> <I>(double tau, const gsl_vector * v, gsl_matrix * A)</I>
<DD><A NAME="IDX1253"></A>
This function applies the Householder matrix P defined by the
scalar <VAR>tau</VAR> and the vector <VAR>v</VAR> to the left-hand side of the
matrix <VAR>A</VAR>. On output the result P A is stored in <VAR>A</VAR>.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_householder_mh</B> <I>(double tau, const gsl_vector * v, gsl_matrix * A)</I>
<DD><A NAME="IDX1254"></A>
This function applies the Householder matrix P defined by the
scalar <VAR>tau</VAR> and the vector <VAR>v</VAR> to the right-hand side of the
matrix <VAR>A</VAR>. On output the result A P is stored in <VAR>A</VAR>.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_householder_hv</B> <I>(double tau, const gsl_vector * v, gsl_vector * w)</I>
<DD><A NAME="IDX1255"></A>
This function applies the Householder transformation P defined by
the scalar <VAR>tau</VAR> and the vector <VAR>v</VAR> to the vector <VAR>w</VAR>. On
output the result P w is stored in <VAR>w</VAR>.
</DL>
</P>
<H2><A NAME="SEC229" HREF="gsl-ref_toc.html#TOC229">Householder solver for linear systems</A></H2>
<P>
<A NAME="IDX1256"></A>
<A NAME="IDX1257"></A>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_HH_solve</B> <I>(gsl_matrix * <VAR>A</VAR>, const gsl_vector * <VAR>b</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1258"></A>
This function solves the system A x = b directly using
Householder transformations. On output the solution is stored in <VAR>x</VAR>
and <VAR>b</VAR> is not modified. The matrix <VAR>A</VAR> is destroyed by the
Householder transformations.
</DL>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_HH_svx</B> <I>(gsl_matrix * <VAR>A</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1259"></A>
This function solves the system A x = b in-place using
Householder transformations. On input <VAR>x</VAR> should contain the
right-hand side b, which is replaced by the solution on output. The
matrix <VAR>A</VAR> is destroyed by the Householder transformations.
</DL>
</P>
<H2><A NAME="SEC230" HREF="gsl-ref_toc.html#TOC230">Tridiagonal Systems</A></H2>
<P>
<A NAME="IDX1260"></A>
</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_solve_tridiag</B> <I>(const gsl_vector * <VAR>diag</VAR>, const gsl_vector * <VAR>e</VAR>, const gsl_vector * <VAR>f</VAR>, const gsl_vector * <VAR>b</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1261"></A>
This function solves the general N-by-N system A x =
b where <VAR>A</VAR> is tridiagonal (
N >= 2). The super-diagonal and
sub-diagonal vectors <VAR>e</VAR> and <VAR>f</VAR> must be one element shorter
than the diagonal vector <VAR>diag</VAR>. The form of <VAR>A</VAR> for the 4-by-4
case is shown below,
</P>
<PRE class="example">
A = ( d_0 e_0 0 0 )
( f_0 d_1 e_1 0 )
( 0 f_1 d_2 e_2 )
( 0 0 f_2 d_3 )
</PRE>
</DL>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_solve_symm_tridiag</B> <I>(const gsl_vector * <VAR>diag</VAR>, const gsl_vector * <VAR>e</VAR>, const gsl_vector * <VAR>b</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1262"></A>
This function solves the general N-by-N system A x =
b where <VAR>A</VAR> is symmetric tridiagonal (
N >= 2). The off-diagonal vector
<VAR>e</VAR> must be one element shorter than the diagonal vector <VAR>diag</VAR>.
The form of <VAR>A</VAR> for the 4-by-4 case is shown below,
</P>
<PRE class="example">
A = ( d_0 e_0 0 0 )
( e_0 d_1 e_1 0 )
( 0 e_1 d_2 e_2 )
( 0 0 e_2 d_3 )
</PRE>
</DL>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_solve_cyc_tridiag</B> <I>(const gsl_vector * <VAR>diag</VAR>, const gsl_vector * <VAR>e</VAR>, const gsl_vector * <VAR>f</VAR>, const gsl_vector * <VAR>b</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1263"></A>
This function solves the general N-by-N system A x =
b where <VAR>A</VAR> is cyclic tridiagonal (
N >= 3). The cyclic super-diagonal and
sub-diagonal vectors <VAR>e</VAR> and <VAR>f</VAR> must have the same number of
elements as the diagonal vector <VAR>diag</VAR>. The form of <VAR>A</VAR> for the
4-by-4 case is shown below,
</P>
<PRE class="example">
A = ( d_0 e_0 0 f_3 )
( f_0 d_1 e_1 0 )
( 0 f_1 d_2 e_2 )
( e_3 0 f_2 d_3 )
</PRE>
</DL>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_linalg_solve_symm_cyc_tridiag</B> <I>(const gsl_vector * <VAR>diag</VAR>, const gsl_vector * <VAR>e</VAR>, const gsl_vector * <VAR>b</VAR>, gsl_vector * <VAR>x</VAR>)</I>
<DD><A NAME="IDX1264"></A>
This function solves the general N-by-N system A x =
b where <VAR>A</VAR> is symmetric cyclic tridiagonal (
N >= 3). The cyclic
off-diagonal vector <VAR>e</VAR> must have the same number of elements as the
diagonal vector <VAR>diag</VAR>. The form of <VAR>A</VAR> for the 4-by-4 case is
shown below,
</P>
<PRE class="example">
A = ( d_0 e_0 0 e_3 )
( e_0 d_1 e_1 0 )
( 0 e_1 d_2 e_2 )
( e_3 0 e_2 d_3 )
</PRE>
</DL>
<H2><A NAME="SEC231" HREF="gsl-ref_toc.html#TOC231">Examples</A></H2>
<P>
The following program solves the linear system A x = b. The
system to be solved is,
</P>
<PRE class="example">
[ 0.18 0.60 0.57 0.96 ] [x0] [1.0]
[ 0.41 0.24 0.99 0.58 ] [x1] = [2.0]
[ 0.14 0.30 0.97 0.66 ] [x2] [3.0]
[ 0.51 0.13 0.19 0.85 ] [x3] [4.0]
</PRE>
<P>
and the solution is found using LU decomposition of the matrix A.
</P>
<PRE class="example">
#include <stdio.h>
#include <gsl/gsl_linalg.h>
int
main (void)
{
double a_data[] = { 0.18, 0.60, 0.57, 0.96,
0.41, 0.24, 0.99, 0.58,
0.14, 0.30, 0.97, 0.66,
0.51, 0.13, 0.19, 0.85 };
double b_data[] = { 1.0, 2.0, 3.0, 4.0 };
gsl_matrix_view m
= gsl_matrix_view_array (a_data, 4, 4);
gsl_vector_view b
= gsl_vector_view_array (b_data, 4);
gsl_vector *x = gsl_vector_alloc (4);
int s;
gsl_permutation * p = gsl_permutation_alloc (4);
gsl_linalg_LU_decomp (&m.matrix, p, &s);
gsl_linalg_LU_solve (&m.matrix, p, &b.vector, x);
printf ("x = \n");
gsl_vector_fprintf (stdout, x, "%g");
gsl_permutation_free (p);
return 0;
}
</PRE>
<P>
Here is the output from the program,
</P>
<PRE class="example">
x = -4.05205
-12.6056
1.66091
8.69377
</PRE>
<P>
This can be verified by multiplying the solution x by the
original matrix A using GNU OCTAVE,
<PRE class="example">
octave> A = [ 0.18, 0.60, 0.57, 0.96;
0.41, 0.24, 0.99, 0.58;
0.14, 0.30, 0.97, 0.66;
0.51, 0.13, 0.19, 0.85 ];
octave> x = [ -4.05205; -12.6056; 1.66091; 8.69377];
octave> A * x
ans =
1.0000
2.0000
3.0000
4.0000
</PRE>
<P>
This reproduces the original right-hand side vector, b, in
accordance with the equation A x = b.
</P>
<H2><A NAME="SEC232" HREF="gsl-ref_toc.html#TOC232">References and Further Reading</A></H2>
<P>
Further information on the algorithms described in this section can be
found in the following book,
</P>
<UL class="itemize">
<LI>
G. H. Golub, C. F. Van Loan, <CITE>Matrix Computations</CITE> (3rd Ed, 1996),
Johns Hopkins University Press, ISBN 0-8018-5414-8.
</UL>
<P>
The LAPACK library is described in the following manual,
</P>
<UL class="itemize">
<LI>
<CITE>LAPACK Users' Guide</CITE> (Third Edition, 1999), Published by SIAM,
ISBN 0-89871-447-8.
<A HREF="http://www.netlib.org/lapack">http://www.netlib.org/lapack</A>
</UL>
<P>
The LAPACK source code can be found at the website above, along
with an online copy of the users guide.
</P>
<P>
The Modified Golub-Reinsch algorithm is described in the following paper,
</P>
<UL class="itemize">
<LI>
T.F. Chan, "An Improved Algorithm for Computing the Singular Value
Decomposition", <CITE>ACM Transactions on Mathematical Software</CITE>, 8
(1982), pp 72--83.
</UL>
<P>
The Jacobi algorithm for singular value decomposition is described in
the following papers,
</P>
<UL class="itemize">
<LI>
J.C.Nash, "A one-sided transformation method for the singular value
decomposition and algebraic eigenproblem", <CITE>Computer Journal</CITE>,
Volume 18, Number 1 (1973), p 74--76
<LI>
James Demmel, Kresimir Veselic, "Jacobi's Method is more accurate than
QR", <CITE>Lapack Working Note 15</CITE> (LAWN-15), October 1989. Available
from netlib, <A HREF="http://www.netlib.org/lapack/">http://www.netlib.org/lapack/</A> in the <CODE>lawns</CODE> or
<CODE>lawnspdf</CODE> directories.
</UL>
<P><HR><P>
<p>Go to the <A HREF="gsl-ref_1.html">first</A>, <A HREF="gsl-ref_12.html">previous</A>, <A HREF="gsl-ref_14.html">next</A>, <A HREF="gsl-ref_50.html">last</A> section, <A HREF="gsl-ref_toc.html">table of contents</A>.
</BODY>
</HTML>
|