File: linalg.texi

package info (click to toggle)
gsl-doc 2.3-1
  • links: PTS
  • area: non-free
  • in suites: buster
  • size: 27,748 kB
  • ctags: 15,177
  • sloc: ansic: 235,014; sh: 11,585; makefile: 925
file content (1767 lines) | stat: -rw-r--r-- 75,665 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
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
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645
1646
1647
1648
1649
1650
1651
1652
1653
1654
1655
1656
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
1672
1673
1674
1675
1676
1677
1678
1679
1680
1681
1682
1683
1684
1685
1686
1687
1688
1689
1690
1691
1692
1693
1694
1695
1696
1697
1698
1699
1700
1701
1702
1703
1704
1705
1706
1707
1708
1709
1710
1711
1712
1713
1714
1715
1716
1717
1718
1719
1720
1721
1722
1723
1724
1725
1726
1727
1728
1729
1730
1731
1732
1733
1734
1735
1736
1737
1738
1739
1740
1741
1742
1743
1744
1745
1746
1747
1748
1749
1750
1751
1752
1753
1754
1755
1756
1757
1758
1759
1760
1761
1762
1763
1764
1765
1766
1767
@cindex linear algebra
@cindex solution of linear systems, Ax=b
@cindex matrix factorization
@cindex factorization of matrices

This chapter describes functions for solving linear systems.  The
library provides linear algebra operations which operate directly on
the @code{gsl_vector} and @code{gsl_matrix} objects.  These routines
use the standard algorithms from Golub & Van Loan's @cite{Matrix
Computations} with Level-1 and Level-2 BLAS calls for efficiency.

The functions described in this chapter are declared in the header file
@file{gsl_linalg.h}.


@menu
* LU Decomposition::            
* QR Decomposition::            
* QR Decomposition with Column Pivoting::
* Complete Orthogonal Decomposition::
* Singular Value Decomposition::  
* Cholesky Decomposition::      
* Pivoted Cholesky Decomposition::      
* Modified Cholesky Decomposition::      
* Tridiagonal Decomposition of Real Symmetric Matrices::  
* Tridiagonal Decomposition of Hermitian Matrices::  
* Hessenberg Decomposition of Real Matrices::
* Hessenberg-Triangular Decomposition of Real Matrices::
* Bidiagonalization::           
* Givens Rotations::  
* Householder Transformations::  
* Householder solver for linear systems::  
* Tridiagonal Systems::
* Triangular Systems::
* Balancing::
* Linear Algebra Examples::     
* Linear Algebra References and Further Reading::  
@end menu

@node LU Decomposition
@section LU Decomposition
@cindex LU decomposition

A general @math{N}-by-@math{N} square matrix @math{A} has an @math{LU} decomposition into
upper and lower triangular matrices,
@tex
\beforedisplay
$$
P A = L U
$$
\afterdisplay
@end tex
@ifinfo

@example
P A = L U
@end example

@end ifinfo
@noindent
where @math{P} is a permutation matrix, @math{L} is unit lower
triangular matrix and @math{U} is upper triangular matrix. For square
matrices this decomposition can be used to convert the linear system
@math{A x = b} into a pair of triangular systems (@math{L y = P b},
@math{U x = y}), which can be solved by forward and back-substitution.
Note that the @math{LU} decomposition is valid for singular matrices.

@deftypefun int gsl_linalg_LU_decomp (gsl_matrix * @var{A}, gsl_permutation * @var{p}, int * @var{signum})
@deftypefunx int gsl_linalg_complex_LU_decomp (gsl_matrix_complex * @var{A}, gsl_permutation * @var{p}, int * @var{signum})
These functions factorize the square matrix @var{A} into the @math{LU}
decomposition @math{PA = LU}.  On output the diagonal and upper
triangular part of the input matrix @var{A} contain the matrix
@math{U}. The lower triangular part of the input matrix (excluding the
diagonal) contains @math{L}.  The diagonal elements of @math{L} are
unity, and are not stored.

The permutation matrix @math{P} is encoded in the permutation
@var{p} on output. The @math{j}-th column of the matrix @math{P}
is given by the @math{k}-th column of the identity matrix, where
@math{k = p_j} the
@math{j}-th element of the permutation vector. The sign of the
permutation is given by @var{signum}. It has the value @math{(-1)^n},
where @math{n} is the number of interchanges in the permutation.

The algorithm used in the decomposition is Gaussian Elimination with
partial pivoting (Golub & Van Loan, @cite{Matrix Computations},
Algorithm 3.4.1).
@end deftypefun

@cindex linear systems, solution of
@deftypefun int gsl_linalg_LU_solve (const gsl_matrix * @var{LU}, const gsl_permutation * @var{p}, const gsl_vector * @var{b}, gsl_vector * @var{x})
@deftypefunx int gsl_linalg_complex_LU_solve (const gsl_matrix_complex * @var{LU}, const gsl_permutation * @var{p}, const gsl_vector_complex * @var{b}, gsl_vector_complex * @var{x})
These functions solve the square system @math{A x = b} using the @math{LU}
decomposition of @math{A} into (@var{LU}, @var{p}) given by
@code{gsl_linalg_LU_decomp} or @code{gsl_linalg_complex_LU_decomp} as input.
@end deftypefun

@deftypefun int gsl_linalg_LU_svx (const gsl_matrix * @var{LU}, const gsl_permutation * @var{p}, gsl_vector * @var{x})
@deftypefunx int gsl_linalg_complex_LU_svx (const gsl_matrix_complex * @var{LU}, const gsl_permutation * @var{p}, gsl_vector_complex * @var{x})
These functions solve the square system @math{A x = b} in-place using the
precomputed @math{LU} decomposition of @math{A} into (@var{LU},@var{p}). On input
@var{x} should contain the right-hand side @math{b}, which is replaced
by the solution on output.
@end deftypefun

@cindex refinement of solutions in linear systems
@cindex iterative refinement of solutions in linear systems
@cindex linear systems, refinement of solutions
@deftypefun int gsl_linalg_LU_refine (const gsl_matrix * @var{A}, const gsl_matrix * @var{LU}, const gsl_permutation * @var{p}, const gsl_vector * @var{b}, gsl_vector * @var{x}, gsl_vector * @var{work})
@deftypefunx int gsl_linalg_complex_LU_refine (const gsl_matrix_complex * @var{A}, const gsl_matrix_complex * @var{LU}, const gsl_permutation * @var{p}, const gsl_vector_complex * @var{b}, gsl_vector_complex * @var{x}, gsl_vector_complex * @var{work})
These functions apply an iterative improvement to @var{x}, the solution
of @math{A x = b}, from the precomputed @math{LU} decomposition of @math{A} into
(@var{LU},@var{p}). Additional workspace of length @var{N} is required in @var{work}.
@end deftypefun

@cindex inverse of a matrix, by LU decomposition
@cindex matrix inverse
@deftypefun int gsl_linalg_LU_invert (const gsl_matrix * @var{LU}, const gsl_permutation * @var{p}, gsl_matrix * @var{inverse})
@deftypefunx int gsl_linalg_complex_LU_invert (const gsl_matrix_complex * @var{LU}, const gsl_permutation * @var{p}, gsl_matrix_complex * @var{inverse})
These functions compute the inverse of a matrix @math{A} from its
@math{LU} decomposition (@var{LU},@var{p}), storing the result in the
matrix @var{inverse}. The inverse is computed by solving the system
@math{A x = b} for each column of the identity matrix.  It is preferable
to avoid direct use of the inverse whenever possible, as the linear
solver functions can obtain the same result more efficiently and
reliably (consult any introductory textbook on numerical linear algebra
for details).
@end deftypefun

@cindex determinant of a matrix, by LU decomposition
@cindex matrix determinant
@deftypefun double gsl_linalg_LU_det (gsl_matrix * @var{LU}, int @var{signum})
@deftypefunx gsl_complex gsl_linalg_complex_LU_det (gsl_matrix_complex * @var{LU}, int @var{signum})
These functions compute the determinant of a matrix @math{A} from its
@math{LU} decomposition, @var{LU}. The determinant is computed as the
product of the diagonal elements of @math{U} and the sign of the row
permutation @var{signum}.
@end deftypefun

@cindex logarithm of the determinant of a matrix
@deftypefun double gsl_linalg_LU_lndet (gsl_matrix * @var{LU})
@deftypefunx double gsl_linalg_complex_LU_lndet (gsl_matrix_complex * @var{LU})
These functions compute the logarithm of the absolute value of the
determinant of a matrix @math{A}, @math{\ln|\det(A)|}, from its @math{LU}
decomposition, @var{LU}.  This function may be useful if the direct
computation of the determinant would overflow or underflow.
@end deftypefun

@cindex sign of the determinant of a matrix
@deftypefun int gsl_linalg_LU_sgndet (gsl_matrix * @var{LU}, int @var{signum})
@deftypefunx gsl_complex gsl_linalg_complex_LU_sgndet (gsl_matrix_complex * @var{LU}, int @var{signum})
These functions compute the sign or phase factor of the determinant of a
matrix @math{A}, @math{\det(A)/|\det(A)|}, from its @math{LU} decomposition,
@var{LU}.
@end deftypefun

@node QR Decomposition
@section QR Decomposition
@cindex QR decomposition

A general rectangular @math{M}-by-@math{N} matrix @math{A} has a
@math{QR} decomposition into the product of an orthogonal
@math{M}-by-@math{M} square matrix @math{Q} (where @math{Q^T Q = I}) and
an @math{M}-by-@math{N} right-triangular matrix @math{R},
@tex
\beforedisplay
$$
A = Q R
$$
\afterdisplay
@end tex
@ifinfo

@example
A = Q R
@end example

@end ifinfo
@noindent
This decomposition can be used to convert the linear system @math{A x =
b} into the triangular system @math{R x = Q^T b}, which can be solved by
back-substitution. Another use of the @math{QR} decomposition is to
compute an orthonormal basis for a set of vectors. The first @math{N}
columns of @math{Q} form an orthonormal basis for the range of @math{A},
@math{ran(A)}, when @math{A} has full column rank.

@deftypefun int gsl_linalg_QR_decomp (gsl_matrix * @var{A}, gsl_vector * @var{tau})
This function factorizes the @math{M}-by-@math{N} matrix @var{A} into
the @math{QR} decomposition @math{A = Q R}.  On output the diagonal and
upper triangular part of the input matrix contain the matrix
@math{R}. The vector @var{tau} and the columns of the lower triangular
part of the matrix @var{A} contain the Householder coefficients and
Householder vectors which encode the orthogonal matrix @var{Q}.  The
vector @var{tau} must be of length @math{k=\min(M,N)}. The matrix
@math{Q} is related to these components by, @math{Q = Q_k ... Q_2 Q_1}
where @math{Q_i = I - \tau_i v_i v_i^T} and @math{v_i} is the
Householder vector @math{v_i =
(0,...,1,A(i+1,i),A(i+2,i),...,A(m,i))}. This is the same storage scheme
as used by @sc{lapack}.

The algorithm used to perform the decomposition is Householder QR (Golub
& Van Loan, @cite{Matrix Computations}, Algorithm 5.2.1).
@end deftypefun

@deftypefun int gsl_linalg_QR_solve (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, const gsl_vector * @var{b}, gsl_vector * @var{x})
This function solves the square system @math{A x = b} using the @math{QR}
decomposition of @math{A} held in (@var{QR}, @var{tau}) which must 
have been computed previously with @code{gsl_linalg_QR_decomp}. 
The least-squares solution for 
rectangular systems can be found using @code{gsl_linalg_QR_lssolve}.
@end deftypefun

@deftypefun int gsl_linalg_QR_svx (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, gsl_vector * @var{x})
This function solves the square system @math{A x = b} in-place using
the @math{QR} decomposition of @math{A} held in (@var{QR},@var{tau})
which must have been computed previously by
@code{gsl_linalg_QR_decomp}.  On input @var{x} should contain the
right-hand side @math{b}, which is replaced by the solution on output.
@end deftypefun

@deftypefun int gsl_linalg_QR_lssolve (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, const gsl_vector * @var{b}, gsl_vector * @var{x}, gsl_vector * @var{residual})
This function finds the least squares solution to the overdetermined
system @math{A x = b} where the matrix @var{A} has more rows than
columns.  The least squares solution minimizes the Euclidean norm of the
residual, @math{||Ax - b||}.The routine requires as input 
the @math{QR} decomposition
of @math{A} into (@var{QR}, @var{tau}) given by
@code{gsl_linalg_QR_decomp}.  The solution is returned in @var{x}.  The
residual is computed as a by-product and stored in @var{residual}.
@end deftypefun

@deftypefun int gsl_linalg_QR_QTvec (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, gsl_vector * @var{v})
This function applies the matrix @math{Q^T} encoded in the decomposition
(@var{QR},@var{tau}) to the vector @var{v}, storing the result @math{Q^T
v} in @var{v}.  The matrix multiplication is carried out directly using
the encoding of the Householder vectors without needing to form the full
matrix @math{Q^T}.
@end deftypefun

@deftypefun int gsl_linalg_QR_Qvec (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, gsl_vector * @var{v})
This function applies the matrix @math{Q} encoded in the decomposition
(@var{QR},@var{tau}) to the vector @var{v}, storing the result @math{Q
v} in @var{v}.  The matrix multiplication is carried out directly using
the encoding of the Householder vectors without needing to form the full
matrix @math{Q}.
@end deftypefun

@deftypefun int gsl_linalg_QR_QTmat (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, gsl_matrix * @var{A})
This function applies the matrix @math{Q^T} encoded in the decomposition
(@var{QR},@var{tau}) to the matrix @var{A}, storing the result @math{Q^T
A} in @var{A}.  The matrix multiplication is carried out directly using
the encoding of the Householder vectors without needing to form the full
matrix @math{Q^T}.
@end deftypefun

@deftypefun int gsl_linalg_QR_Rsolve (const gsl_matrix * @var{QR}, const gsl_vector * @var{b}, gsl_vector * @var{x})
This function solves the triangular system @math{R x = b} for
@var{x}. It may be useful if the product @math{b' = Q^T b} has already
been computed using @code{gsl_linalg_QR_QTvec}.
@end deftypefun

@deftypefun int gsl_linalg_QR_Rsvx (const gsl_matrix * @var{QR}, gsl_vector * @var{x})
This function solves the triangular system @math{R x = b} for @var{x}
in-place. On input @var{x} should contain the right-hand side @math{b}
and is replaced by the solution on output. This function may be useful if
the product @math{b' = Q^T b} has already been computed using
@code{gsl_linalg_QR_QTvec}.
@end deftypefun

@deftypefun int gsl_linalg_QR_unpack (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, gsl_matrix * @var{Q}, gsl_matrix * @var{R})
This function unpacks the encoded @math{QR} decomposition
(@var{QR},@var{tau}) into the matrices @var{Q} and @var{R}, where
@var{Q} is @math{M}-by-@math{M} and @var{R} is @math{M}-by-@math{N}. 
@end deftypefun

@deftypefun int gsl_linalg_QR_QRsolve (gsl_matrix * @var{Q}, gsl_matrix * @var{R}, const gsl_vector * @var{b}, gsl_vector * @var{x})
This function solves the system @math{R x = Q^T b} for @var{x}. It can
be used when the @math{QR} decomposition of a matrix is available in
unpacked form as (@var{Q}, @var{R}).
@end deftypefun

@deftypefun int gsl_linalg_QR_update (gsl_matrix * @var{Q}, gsl_matrix * @var{R}, gsl_vector * @var{w}, const gsl_vector * @var{v})
This function performs a rank-1 update @math{w v^T} of the @math{QR}
decomposition (@var{Q}, @var{R}). The update is given by @math{Q'R' = Q
(R + w v^T)} where the output matrices @math{Q'} and @math{R'} are also
orthogonal and right triangular. Note that @var{w} is destroyed by the
update.
@end deftypefun

@deftypefun int gsl_linalg_R_solve (const gsl_matrix * @var{R}, const gsl_vector * @var{b}, gsl_vector * @var{x})
This function solves the triangular system @math{R x = b} for the
@math{N}-by-@math{N} matrix @var{R}.
@end deftypefun

@deftypefun int gsl_linalg_R_svx (const gsl_matrix * @var{R}, gsl_vector * @var{x})
This function solves the triangular system @math{R x = b} in-place. On
input @var{x} should contain the right-hand side @math{b}, which is
replaced by the solution on output.
@end deftypefun

@node QR Decomposition with Column Pivoting
@section QR Decomposition with Column Pivoting
@cindex QR decomposition with column pivoting

The @math{QR} decomposition of an @math{M}-by-@math{N} matrix @math{A}
can be extended to the rank deficient case by introducing a column permutation @math{P},
@tex
\beforedisplay
$$
A P = Q R
$$
\afterdisplay
@end tex
@ifinfo

@example
A P = Q R
@end example

@end ifinfo
@noindent
The first @math{r} columns of @math{Q} form an orthonormal basis
for the range of @math{A} for a matrix with column rank @math{r}.  This
decomposition can also be used to convert the linear system @math{A x =
b} into the triangular system @math{R y = Q^T b, x = P y}, which can be
solved by back-substitution and permutation.  We denote the @math{QR}
decomposition with column pivoting by @math{QRP^T} since @math{A = Q R
P^T}. When @math{A} is rank deficient with @math{r = {\rm rank}(A)}, the matrix
@math{R} can be partitioned as
@tex
\beforedisplay
$$
R = \left(
\matrix{
R_{11} & R_{12} \cr
0 & R_{22} \cr
}
\right) \approx
\left(
\matrix{
R_{11} & R_{12} \cr
0 & 0 \cr
}
\right)
$$
\afterdisplay
@end tex
@ifinfo

@example
R = [ R11 R12; 0 R22 ] =~ [ R11 R12; 0 0 ]
@end example

@end ifinfo
where @math{R_{11}} is @math{r}-by-@math{r} and nonsingular. In this case,
a ``basic'' least squares solution for the overdetermined system @math{A x = b}
can be obtained as
@tex
\beforedisplay
$$
x = P \left(
\matrix{
R_{11}^{-1} c_1 \cr
0 \cr
}
\right)
$$
\afterdisplay
@end tex
@ifinfo

@example
x = P [ R11^-1 c1 ; 0 ]
@end example

@end ifinfo
where @math{c_1} consists of the first @math{r} elements of @math{Q^T b}.
This basic solution is not guaranteed to be the minimum norm solution unless
@math{R_{12} = 0} (see @ref{Complete Orthogonal Decomposition}).

@deftypefun int gsl_linalg_QRPT_decomp (gsl_matrix * @var{A}, gsl_vector * @var{tau}, gsl_permutation * @var{p}, int * @var{signum}, gsl_vector * @var{norm})
This function factorizes the @math{M}-by-@math{N} matrix @var{A} into
the @math{QRP^T} decomposition @math{A = Q R P^T}.  On output the
diagonal and upper triangular part of the input matrix contain the
matrix @math{R}. The permutation matrix @math{P} is stored in the
permutation @var{p}.  The sign of the permutation is given by
@var{signum}. It has the value @math{(-1)^n}, where @math{n} is the
number of interchanges in the permutation. The vector @var{tau} and the
columns of the lower triangular part of the matrix @var{A} contain the
Householder coefficients and vectors which encode the orthogonal matrix
@var{Q}.  The vector @var{tau} must be of length @math{k=\min(M,N)}. The
matrix @math{Q} is related to these components by, @math{Q = Q_k ... Q_2
Q_1} where @math{Q_i = I - \tau_i v_i v_i^T} and @math{v_i} is the
Householder vector @math{v_i =
(0,...,1,A(i+1,i),A(i+2,i),...,A(m,i))}. This is the same storage scheme
as used by @sc{lapack}.  The vector @var{norm} is a workspace of length
@var{N} used for column pivoting.

The algorithm used to perform the decomposition is Householder QR with
column pivoting (Golub & Van Loan, @cite{Matrix Computations}, Algorithm
5.4.1).
@end deftypefun

@deftypefun int gsl_linalg_QRPT_decomp2 (const gsl_matrix * @var{A}, gsl_matrix * @var{q}, gsl_matrix * @var{r}, gsl_vector * @var{tau}, gsl_permutation * @var{p}, int * @var{signum}, gsl_vector * @var{norm})
This function factorizes the matrix @var{A} into the decomposition
@math{A = Q R P^T} without modifying @var{A} itself and storing the
output in the separate matrices @var{q} and @var{r}.
@end deftypefun

@deftypefun int gsl_linalg_QRPT_solve (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, const gsl_permutation * @var{p}, const gsl_vector * @var{b}, gsl_vector * @var{x})
This function solves the square system @math{A x = b} using the @math{QRP^T}
decomposition of @math{A} held in (@var{QR}, @var{tau}, @var{p}) which must 
have been computed previously by @code{gsl_linalg_QRPT_decomp}.
@end deftypefun

@deftypefun int gsl_linalg_QRPT_svx (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, const gsl_permutation * @var{p}, gsl_vector * @var{x})
This function solves the square system @math{A x = b} in-place using the
@math{QRP^T} decomposition of @math{A} held in
(@var{QR},@var{tau},@var{p}). On input @var{x} should contain the
right-hand side @math{b}, which is replaced by the solution on output.
@end deftypefun

@deftypefun int gsl_linalg_QRPT_lssolve (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, const gsl_permutation * @var{p}, const gsl_vector * @var{b}, gsl_vector * @var{x}, gsl_vector * @var{residual})
This function finds the least squares solution to the overdetermined
system @math{A x = b} where the matrix @var{A} has more rows than
columns and is assumed to have full rank. The least squares solution minimizes
the Euclidean norm of the residual, @math{||b - A x||}. The routine requires as input 
the @math{QR} decomposition of @math{A} into (@var{QR}, @var{tau}, @var{p}) given by
@code{gsl_linalg_QRPT_decomp}.  The solution is returned in @var{x}.  The
residual is computed as a by-product and stored in @var{residual}. For rank
deficient matrices, @code{gsl_linalg_QRPT_lssolve2} should be used instead.
@end deftypefun

@deftypefun int gsl_linalg_QRPT_lssolve2 (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, const gsl_permutation * @var{p}, const gsl_vector * @var{b}, const size_t @var{rank}, gsl_vector * @var{x}, gsl_vector * @var{residual})
This function finds the least squares solution to the overdetermined
system @math{A x = b} where the matrix @var{A} has more rows than
columns and has rank given by the input @var{rank}. If the user does not
know the rank of @math{A}, the routine @code{gsl_linalg_QRPT_rank} can be
called to estimate it. The least squares solution is
the so-called ``basic'' solution discussed above and may not be the minimum
norm solution. The routine requires as input 
the @math{QR} decomposition of @math{A} into (@var{QR}, @var{tau}, @var{p}) given by
@code{gsl_linalg_QRPT_decomp}.  The solution is returned in @var{x}.  The
residual is computed as a by-product and stored in @var{residual}.
@end deftypefun

@deftypefun int gsl_linalg_QRPT_QRsolve (const gsl_matrix * @var{Q}, const gsl_matrix * @var{R}, const gsl_permutation * @var{p}, const gsl_vector * @var{b}, gsl_vector * @var{x})
This function solves the square system @math{R P^T x = Q^T b} for
@var{x}. It can be used when the @math{QR} decomposition of a matrix is
available in unpacked form as (@var{Q}, @var{R}).
@end deftypefun

@deftypefun int gsl_linalg_QRPT_update (gsl_matrix * @var{Q}, gsl_matrix * @var{R}, const gsl_permutation * @var{p}, gsl_vector * @var{w}, const gsl_vector * @var{v})
This function performs a rank-1 update @math{w v^T} of the @math{QRP^T}
decomposition (@var{Q}, @var{R}, @var{p}). The update is given by
@math{Q'R' = Q (R + w v^T P)} where the output matrices @math{Q'} and
@math{R'} are also orthogonal and right triangular. Note that @var{w} is
destroyed by the update. The permutation @var{p} is not changed.
@end deftypefun

@deftypefun int gsl_linalg_QRPT_Rsolve (const gsl_matrix * @var{QR}, const gsl_permutation * @var{p}, const gsl_vector * @var{b}, gsl_vector * @var{x})
This function solves the triangular system @math{R P^T x = b} for the
@math{N}-by-@math{N} matrix @math{R} contained in @var{QR}.
@end deftypefun

@deftypefun int gsl_linalg_QRPT_Rsvx (const gsl_matrix * @var{QR}, const gsl_permutation * @var{p}, gsl_vector * @var{x})
This function solves the triangular system @math{R P^T x = b} in-place
for the @math{N}-by-@math{N} matrix @math{R} contained in @var{QR}. On
input @var{x} should contain the right-hand side @math{b}, which is
replaced by the solution on output.
@end deftypefun

@deftypefun size_t gsl_linalg_QRPT_rank (const gsl_matrix * @var{QR}, const double @var{tol})
This function estimates the rank of the triangular matrix @math{R} contained in @var{QR}.
The algorithm simply counts the number of diagonal elements of @math{R} whose absolute value
is greater than the specified tolerance @var{tol}. If the input @var{tol} is negative,
a default value of @math{20 (M + N) eps(max(|diag(R)|))} is used.
@end deftypefun

@deftypefun int gsl_linalg_QRPT_rcond (const gsl_matrix * @var{QR}, double * @var{rcond}, gsl_vector * @var{work})
This function estimates the reciprocal condition number (using the 1-norm) of the @math{R} factor,
stored in the upper triangle of @var{QR}. The reciprocal condition number estimate, defined as
@math{1 / (||R||_1 \cdot ||R^{-1}||_1)}, is stored in @var{rcond}.
Additional workspace of size @math{3 N} is required in @var{work}.
@end deftypefun

@node Complete Orthogonal Decomposition
@section Complete Orthogonal Decomposition

The complete orthogonal decomposition of a @math{M}-by-@math{N} matrix
@math{A} is a generalization of the QR decomposition with column pivoting, given by
@tex
\beforedisplay
$$
A P = Q
\left(
\matrix{
R_{11} & 0 \cr
0 & 0 \cr
}
\right) Z
$$
\afterdisplay
@end tex
@ifinfo

@example
A P = Q [ R11 0 ] Z
        [  0  0 ]
@end example

@end ifinfo
@noindent
where @math{P} is a @math{N}-by-@math{N} permutation matrix,
@math{Q} is @math{M}-by-@math{M} orthogonal, @math{R_{11}} is
@math{r}-by-@math{r} upper triangular, with @math{r = {\rm rank}(A)},
and @math{Z} is @math{N}-by-@math{N} orthogonal. If @math{A}
has full rank, then @math{R_{11} = R}, @math{Z = I} and this reduces to the
QR decomposition with column pivoting. The advantage of using
the complete orthogonal decomposition for rank deficient matrices
is the ability to compute the minimum norm solution to the linear
least squares problem @math{Ax = b}, which is given by
@tex
\beforedisplay
$$
x = P Z^T
\left(
\matrix{
R_{11}^{-1} c_1 \cr
0 \cr
}
\right)
$$
\afterdisplay
@end tex
@ifinfo

@example
x = P Z^T [ R11^-1 c1 ]
          [    0      ]
@end example

@end ifinfo
and the vector @math{c_1} is the first @math{r} elements of @math{Q^T b}.

@deftypefun int gsl_linalg_COD_decomp (gsl_matrix * @var{A}, gsl_vector * @var{tau_Q}, gsl_vector * @var{tau_Z}, gsl_permutation * @var{p}, size_t * @var{rank}, gsl_vector * @var{work})
@deftypefunx int gsl_linalg_COD_decomp_e (gsl_matrix * @var{A}, gsl_vector * @var{tau_Q}, gsl_vector * @var{tau_Z}, gsl_permutation * @var{p}, double @var{tol}, size_t * @var{rank}, gsl_vector * @var{work})
These functions factor the @math{M}-by-@math{N} matrix @var{A} into the decomposition @math{A = Q R Z P^T}. The rank of @var{A}
is computed as the number of diagonal elements of @math{R} greater than the tolerance @var{tol} and output in @var{rank}.
If @var{tol} is not specified, a default value is used (see @code{gsl_linalg_QRPT_rank}). On output, the permutation
matrix @math{P} is stored in @var{p}. The matrix @math{R_{11}} is stored in the upper @var{rank}-by-@var{rank} block of @var{A}.
The matrices @math{Q} and @math{Z} are encoded in packed storage in @var{A} on output. The vectors @var{tau_Q} and @var{tau_Z}
contain the Householder scalars corresponding to the matrices @math{Q} and @math{Z} respectively and must be
of length @math{k = \min(M,N)}. The vector @var{work} is additional workspace of length @math{N}.
@end deftypefun

@deftypefun int gsl_linalg_COD_lssolve (const gsl_matrix * @var{QRZ}, const gsl_vector * @var{tau_Q}, const gsl_vector * @var{tau_Z}, const gsl_permutation * @var{p}, const size_t @var{rank}, const gsl_vector * @var{b}, gsl_vector * @var{x}, gsl_vector * @var{residual})
This function finds the least squares solution to the overdetermined
system @math{A x = b} where the matrix @var{A} has more rows than
columns.  The least squares solution minimizes the Euclidean norm of the
residual, @math{||b - A x||}.  The routine requires as input 
the @math{QRZ} decomposition of @math{A} into (@var{QRZ}, @var{tau_Q}, @var{tau_Z}, @var{p}, @var{rank})
given by @code{gsl_linalg_COD_decomp}.  The solution is returned in @var{x}.  The
residual is computed as a by-product and stored in @var{residual}.
@end deftypefun

@deftypefun int gsl_linalg_COD_unpack (const gsl_matrix * @var{QRZ}, const gsl_vector * @var{tau_Q}, const gsl_vector * @var{tau_Z}, const size_t @var{rank}, gsl_matrix * @var{Q}, gsl_matrix * @var{R}, gsl_matrix * @var{Z})
This function unpacks the encoded @math{QRZ} decomposition
(@var{QRZ}, @var{tau_Q}, @var{tau_Z}, @var{rank}) into the matrices
@var{Q}, @var{R}, and @var{Z}, where @var{Q} is @math{M}-by-@math{M},
@var{R} is @math{M}-by-@math{N}, and @var{Z} is @math{N}-by-@math{N}.
@end deftypefun

@deftypefun int gsl_linalg_COD_matZ (const gsl_matrix * @var{QRZ}, const gsl_vector * @var{tau_Z}, const size_t @var{rank}, gsl_matrix * @var{A}, gsl_vector * @var{work})
This function multiplies the input matrix @var{A} on the right by @var{Z},
@math{A' = A Z} using the encoded @math{QRZ} decomposition
(@var{QRZ}, @var{tau_Z}, @var{rank}). @var{A} must have @math{N} columns but may
have any number of rows. Additional workspace of length @math{M} is provided
in @var{work}.
@end deftypefun

@node Singular Value Decomposition
@section Singular Value Decomposition
@cindex SVD
@cindex singular value decomposition

A general rectangular @math{M}-by-@math{N} matrix @math{A} has a
singular value decomposition (@sc{svd}) into the product of an
@math{M}-by-@math{N} orthogonal matrix @math{U}, an @math{N}-by-@math{N}
diagonal matrix of singular values @math{S} and the transpose of an
@math{N}-by-@math{N} orthogonal square matrix @math{V},
@tex
\beforedisplay
$$
A = U S V^T
$$
\afterdisplay
@end tex
@ifinfo

@example
A = U S V^T
@end example

@end ifinfo
@noindent
The singular values
@c{$\sigma_i = S_{ii}$}
@math{\sigma_i = S_@{ii@}} are all non-negative and are
generally chosen to form a non-increasing sequence 
@c{$\sigma_1 \ge \sigma_2 \ge ... \ge \sigma_N \ge 0$}
@math{\sigma_1 >= \sigma_2 >= ... >= \sigma_N >= 0}.  

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.

For a rank-deficient matrix, the null space of @math{A} is given by
the columns of @math{V} corresponding to the zero singular values.
Similarly, the range of @math{A} is given by columns of @math{U}
corresponding to the non-zero singular values.

Note that the routines here compute the ``thin'' version of the SVD
with @math{U} as @math{M}-by-@math{N} orthogonal matrix. This allows
in-place computation and is the most commonly-used form in practice.
Mathematically, the ``full'' SVD is defined with @math{U} as an
@math{M}-by-@math{M} orthogonal matrix and @math{S} as an
@math{M}-by-@math{N} diagonal matrix (with additional rows of zeros).

@deftypefun int gsl_linalg_SV_decomp (gsl_matrix * @var{A}, gsl_matrix * @var{V}, gsl_vector * @var{S}, gsl_vector * @var{work})
This function factorizes the @math{M}-by-@math{N} matrix @var{A} into
the singular value decomposition @math{A = U S V^T} for @c{$M \ge N$}
@math{M >= N}.  On output the matrix @var{A} is replaced by
@math{U}. The diagonal elements of the singular value matrix @math{S}
are stored in the vector @var{S}. The singular values are non-negative
and form a non-increasing sequence from @math{S_1} to @math{S_N}. The
matrix @var{V} contains the elements of @math{V} in untransposed
form. To form the product @math{U S V^T} it is necessary to take the
transpose of @var{V}.  A workspace of length @var{N} is required in
@var{work}.

This routine uses the Golub-Reinsch SVD algorithm.  
@end deftypefun

@deftypefun int gsl_linalg_SV_decomp_mod (gsl_matrix * @var{A}, gsl_matrix * @var{X}, gsl_matrix * @var{V}, gsl_vector * @var{S}, gsl_vector * @var{work})
This function computes the SVD using the modified Golub-Reinsch
algorithm, which is faster for @c{$M \gg N$}
@math{M>>N}.  It requires the vector @var{work} of length @var{N} and the
@math{N}-by-@math{N} matrix @var{X} as additional working space.
@end deftypefun

@deftypefun int gsl_linalg_SV_decomp_jacobi (gsl_matrix * @var{A}, gsl_matrix * @var{V}, gsl_vector * @var{S})
@cindex Jacobi orthogonalization
This function computes the SVD of the @math{M}-by-@math{N} matrix @var{A}
using one-sided Jacobi orthogonalization for @c{$M \ge N$} 
@math{M >= N}.  The Jacobi method can compute singular values to higher
relative accuracy than Golub-Reinsch algorithms (see references for
details).
@end deftypefun

@deftypefun int gsl_linalg_SV_solve (const gsl_matrix * @var{U}, const gsl_matrix * @var{V}, const gsl_vector * @var{S}, const gsl_vector * @var{b}, gsl_vector * @var{x})
This function solves the system @math{A x = b} using the singular value
decomposition (@var{U}, @var{S}, @var{V}) of @math{A} which must 
have been computed previously with @code{gsl_linalg_SV_decomp}.

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. 

In the over-determined case where @var{A} has more rows than columns the
system is solved in the least squares sense, returning the solution
@var{x} which minimizes @math{||A x - b||_2}.
@end deftypefun

@deftypefun int gsl_linalg_SV_leverage (const gsl_matrix * @var{U}, gsl_vector * @var{h})
This function computes the statistical leverage values @math{h_i} of a matrix @math{A}
using its singular value decomposition (@var{U}, @var{S}, @var{V}) previously computed
with @code{gsl_linalg_SV_decomp}. @math{h_i} are the diagonal values of the matrix
@math{A (A^T A)^{-1} A^T} and depend only on the matrix @var{U} which is the input to
this function.
@end deftypefun

@node Cholesky Decomposition
@section Cholesky Decomposition
@cindex Cholesky decomposition
@cindex square root of a matrix, Cholesky decomposition
@cindex matrix square root, Cholesky decomposition

A symmetric, positive definite square matrix @math{A} has a Cholesky
decomposition into a product of a lower triangular matrix @math{L} and
its transpose @math{L^T},
@tex
\beforedisplay
$$
A = L L^T
$$
\afterdisplay
@end tex
@ifinfo

@example
A = L L^T
@end example

@end ifinfo
@noindent
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 @math{A x = b} into a pair of triangular systems
(@math{L y = b}, @math{L^T x = y}), which can be solved by forward and
back-substitution.

If the matrix @math{A} is near singular, it is sometimes possible to reduce
the condition number and recover a more accurate solution vector @math{x}
by scaling as
@tex
\beforedisplay
$$
\left( S A S \right) \left( S^{-1} x \right) = S b
$$
\afterdisplay
@end tex
@ifinfo

@example
( S A S ) ( S^(-1) x ) = S b
@end example

@end ifinfo
where @math{S} is a diagonal matrix whose elements are given by
@math{S_{ii} = 1/\sqrt{A_{ii}}}. This scaling is also known as
Jacobi preconditioning. There are routines below to solve
both the scaled and unscaled systems.

@deftypefun int gsl_linalg_cholesky_decomp1 (gsl_matrix * @var{A})
@deftypefunx int gsl_linalg_complex_cholesky_decomp (gsl_matrix_complex * @var{A})
These functions factorize the symmetric, positive-definite square matrix
@var{A} into the Cholesky decomposition @math{A = L L^T} (or
@c{$A = L L^{\dagger}$}
@math{A = L L^H}
for the complex case). On input, the values from the diagonal and lower-triangular
part of the matrix @var{A} are used (the upper triangular part is ignored).  On output
the diagonal and lower triangular part of the input matrix @var{A} contain the matrix
@math{L}, while the upper triangular part is unmodified.  If the matrix is not
positive-definite then the decomposition will fail, returning the
error code @code{GSL_EDOM}.  

When testing whether a matrix is positive-definite, disable the error
handler first to avoid triggering an error.
@end deftypefun

@deftypefun int gsl_linalg_cholesky_decomp (gsl_matrix * @var{A})
This function is now deprecated and is provided only for backward compatibility.
@end deftypefun

@deftypefun int gsl_linalg_cholesky_solve (const gsl_matrix * @var{cholesky}, const gsl_vector * @var{b}, gsl_vector * @var{x})
@deftypefunx int gsl_linalg_complex_cholesky_solve (const gsl_matrix_complex * @var{cholesky}, const gsl_vector_complex * @var{b}, gsl_vector_complex * @var{x})
These functions solve the system @math{A x = b} using the Cholesky
decomposition of @math{A} held in the matrix @var{cholesky} which must
have been previously computed by @code{gsl_linalg_cholesky_decomp} or
@code{gsl_linalg_complex_cholesky_decomp}.
@end deftypefun

@deftypefun int gsl_linalg_cholesky_svx (const gsl_matrix * @var{cholesky}, gsl_vector * @var{x})
@deftypefunx int gsl_linalg_complex_cholesky_svx (const gsl_matrix_complex * @var{cholesky}, gsl_vector_complex * @var{x})
These functions solve the system @math{A x = b} in-place using the
Cholesky decomposition of @math{A} held in the matrix @var{cholesky}
which must have been previously computed by
@code{gsl_linalg_cholesky_decomp} or
@code{gsl_linalg_complex_cholesky_decomp}. On input @var{x} should
contain the right-hand side @math{b}, which is replaced by the
solution on output.
@end deftypefun

@deftypefun int gsl_linalg_cholesky_invert (gsl_matrix * @var{cholesky})
@deftypefunx int gsl_linalg_complex_cholesky_invert (gsl_matrix_complex * @var{cholesky})
These functions compute the inverse of a matrix from its Cholesky
decomposition @var{cholesky}, which must have been previously computed
by @code{gsl_linalg_cholesky_decomp} or
@code{gsl_linalg_complex_cholesky_decomp}.  On output, the inverse is
stored in-place in @var{cholesky}.
@end deftypefun

@deftypefun int gsl_linalg_cholesky_decomp2 (gsl_matrix * @var{A}, gsl_vector * @var{S})
This function calculates a diagonal scaling transformation @math{S} for
the symmetric, positive-definite square matrix @var{A}, and then
computes the Cholesky decomposition @math{S A S = L L^T}.
On input, the values from the diagonal and lower-triangular part of the matrix @var{A} are
used (the upper triangular part is ignored).  On output the diagonal and lower triangular part
of the input matrix @var{A} contain the matrix @math{L}, while the upper triangular part
of the input matrix is overwritten with @math{L^T} (the diagonal terms being
identical for both @math{L} and @math{L^T}).  If the matrix is not
positive-definite then the decomposition will fail, returning the
error code @code{GSL_EDOM}. The diagonal scale factors are stored in @var{S}
on output.

When testing whether a matrix is positive-definite, disable the error
handler first to avoid triggering an error.
@end deftypefun

@deftypefun int gsl_linalg_cholesky_solve2 (const gsl_matrix * @var{cholesky}, const gsl_vector * @var{S}, const gsl_vector * @var{b}, gsl_vector * @var{x})
This function solves the system @math{(S A S) (S^{-1} x) = S b} using the Cholesky
decomposition of @math{S A S} held in the matrix @var{cholesky} which must
have been previously computed by @code{gsl_linalg_cholesky_decomp2}.
@end deftypefun

@deftypefun int gsl_linalg_cholesky_svx2 (const gsl_matrix * @var{cholesky}, const gsl_vector * @var{S}, gsl_vector * @var{x})
This function solves the system @math{(S A S) (S^{-1} x) = S b} in-place using the
Cholesky decomposition of @math{S A S} held in the matrix @var{cholesky}
which must have been previously computed by
@code{gsl_linalg_cholesky_decomp2}.  On input @var{x} should
contain the right-hand side @math{b}, which is replaced by the
solution on output.
@end deftypefun

@deftypefun int gsl_linalg_cholesky_scale (const gsl_matrix * @var{A}, gsl_vector * @var{S})
This function calculates a diagonal scaling transformation of the
symmetric, positive definite matrix @var{A}, such that
@math{S A S} has a condition number within a factor of @math{N}
of the matrix of smallest possible condition number over all
possible diagonal scalings. On output, @var{S} contains the
scale factors, given by @math{S_i = 1/\sqrt{A_{ii}}}.
For any @math{A_{ii} \le 0}, the corresponding scale factor @math{S_i}
is set to @math{1}.
@end deftypefun

@deftypefun int gsl_linalg_cholesky_scale_apply (gsl_matrix * @var{A}, const gsl_vector * @var{S})
This function applies the scaling transformation @var{S} to the matrix @var{A}. On output,
@var{A} is replaced by @math{S A S}.
@end deftypefun

@deftypefun int gsl_linalg_cholesky_rcond (const gsl_matrix * @var{cholesky}, double * @var{rcond}, gsl_vector * @var{work})
This function estimates the reciprocal condition number (using the 1-norm) of the symmetric positive
definite matrix @math{A}, using its Cholesky decomposition provided in @var{cholesky}.
The reciprocal condition number estimate, defined as @math{1 / (||A||_1 \cdot ||A^{-1}||_1)}, is stored
in @var{rcond}.  Additional workspace of size @math{3 N} is required in @var{work}.
@end deftypefun

@node Pivoted Cholesky Decomposition
@section Pivoted Cholesky Decomposition
@cindex Cholesky decomposition, pivoted
@cindex Pivoted Cholesky Decomposition

A symmetric, positive definite square matrix @math{A} has an alternate
Cholesky decomposition into a product of a lower unit triangular matrix @math{L},
a diagonal matrix @math{D} and @math{L^T}, given by @math{L D L^T}. This is
equivalent to the Cholesky formulation discussed above, with
the standard Cholesky lower triangular factor given by @math{L D^{1 \over 2}}.
For ill-conditioned matrices, it can help to use a pivoting strategy to
prevent the entries of @math{D} and @math{L} from growing too large, and also
ensure @math{D_1 \ge D_2 \ge \cdots \ge D_n > 0}, where @math{D_i} are
the diagonal entries of @math{D}. The final decomposition is given by
@tex
\beforedisplay
$$
P A P^T = L D L^T
$$
\afterdisplay
@end tex
@ifinfo

@example
P A P^T = L D L^T
@end example

@end ifinfo
where @math{P} is a permutation matrix.

@deftypefun int gsl_linalg_pcholesky_decomp (gsl_matrix * @var{A}, gsl_permutation * @var{p})
This function factors the symmetric, positive-definite square matrix
@var{A} into the Pivoted Cholesky decomposition @math{P A P^T = L D L^T}.
On input, the values from the diagonal and lower-triangular part of the matrix @var{A} are
used to construct the factorization. On output the diagonal of the input matrix @var{A} stores
the diagonal elements of @math{D}, and the lower triangular portion of @var{A}
contains the matrix @math{L}. Since @math{L} has ones on its diagonal these
do not need to be explicitely stored. The upper triangular portion of @var{A} is
unmodified. The permutation matrix @math{P} is stored in @var{p} on output.
@end deftypefun

@deftypefun int gsl_linalg_pcholesky_solve (const gsl_matrix * @var{LDLT}, const gsl_permutation * @var{p}, const gsl_vector * @var{b}, gsl_vector * @var{x})
This function solves the system @math{A x = b} using the Pivoted Cholesky
decomposition of @math{A} held in the matrix @var{LDLT} and permutation
@var{p} which must have been previously computed by @code{gsl_linalg_pcholesky_decomp}.
@end deftypefun

@deftypefun int gsl_linalg_pcholesky_svx (const gsl_matrix * @var{LDLT}, const gsl_permutation * @var{p}, gsl_vector * @var{x})
This function solves the system @math{A x = b} in-place using the Pivoted Cholesky
decomposition of @math{A} held in the matrix @var{LDLT} and permutation
@var{p} which must have been previously computed by @code{gsl_linalg_pcholesky_decomp}.
On input, @var{x} contains the right hand side vector @math{b} which is
replaced by the solution vector on output.
@end deftypefun

@deftypefun int gsl_linalg_pcholesky_decomp2 (gsl_matrix * @var{A}, gsl_permutation * @var{p}, gsl_vector * @var{S})
This function computes the pivoted Cholesky factorization of the matrix
@math{S A S}, where the input matrix @var{A} is symmetric and positive
definite, and the diagonal scaling matrix @var{S} is computed to reduce the
condition number of @var{A} as much as possible. See
@ref{Cholesky Decomposition} for more information on the matrix @var{S}.
The Pivoted Cholesky decomposition satisfies @math{P S A S P^T = L D L^T}.
On input, the values from the diagonal and lower-triangular part of the matrix @var{A} are
used to construct the factorization.  On output the diagonal of the input matrix @var{A} stores the diagonal
elements of @math{D}, and the lower triangular portion of @var{A}
contains the matrix @math{L}. Since @math{L} has ones on its diagonal these
do not need to be explicitely stored. The upper triangular portion of @var{A}
is unmodified. The permutation matrix @math{P} is stored in @var{p} on output.
The diagonal scaling transformation is stored in @var{S} on output.
@end deftypefun

@deftypefun int gsl_linalg_pcholesky_solve2 (const gsl_matrix * @var{LDLT}, const gsl_permutation * @var{p}, const gsl_vector * @var{S}, const gsl_vector * @var{b}, gsl_vector * @var{x})
This function solves the system @math{(S A S) (S^{-1} x) = S b} using the Pivoted Cholesky
decomposition of @math{S A S} held in the matrix @var{LDLT}, permutation
@var{p}, and vector @var{S}, which must have been previously computed by
@code{gsl_linalg_pcholesky_decomp2}.
@end deftypefun

@deftypefun int gsl_linalg_pcholesky_svx2 (const gsl_matrix * @var{LDLT}, const gsl_permutation * @var{p}, const gsl_vector * @var{S}, gsl_vector * @var{x})
This function solves the system @math{(S A S) (S^{-1} x) = S b} in-place using the Pivoted Cholesky
decomposition of @math{S A S} held in the matrix @var{LDLT}, permutation
@var{p} and vector @var{S}, which must have been previously computed by
@code{gsl_linalg_pcholesky_decomp2}.
On input, @var{x} contains the right hand side vector @math{b} which is
replaced by the solution vector on output.
@end deftypefun

@deftypefun int gsl_linalg_pcholesky_invert (const gsl_matrix * @var{LDLT}, const gsl_permutation * @var{p}, gsl_matrix * @var{Ainv})
This function computes the inverse of the matrix @math{A}, using the Pivoted
Cholesky decomposition stored in @var{LDLT} and @var{p}. On output, the
matrix @var{Ainv} contains @math{A^{-1}}.
@end deftypefun

@deftypefun int gsl_linalg_pcholesky_rcond (const gsl_matrix * @var{LDLT}, const gsl_permutation * @var{p}, double * @var{rcond}, gsl_vector * @var{work})
This function estimates the reciprocal condition number (using the 1-norm) of the symmetric positive
definite matrix @math{A}, using its pivoted Cholesky decomposition provided in @var{LDLT}.
The reciprocal condition number estimate, defined as @math{1 / (||A||_1 \cdot ||A^{-1}||_1)}, is stored
in @var{rcond}.  Additional workspace of size @math{3 N} is required in @var{work}.
@end deftypefun

@node Modified Cholesky Decomposition
@section Modified Cholesky Decomposition
@cindex Cholesky decomposition, modified
@cindex Modified Cholesky Decomposition

The modified Cholesky decomposition is suitable for solving systems
@math{A x = b} where @math{A} is a symmetric indefinite matrix. Such
matrices arise in nonlinear optimization algorithms. The standard
Cholesky decomposition requires a positive definite matrix and would
fail in this case. Instead of resorting to a method like QR or SVD,
which do not take into account the symmetry of the matrix, we can
instead introduce a small perturbation to the matrix @math{A} to
make it positive definite, and then use a Cholesky decomposition on
the perturbed matrix. The resulting decomposition satisfies
@tex
\beforedisplay
$$
P (A + E) P^T = L D L^T
$$
\afterdisplay
@end tex
@ifinfo

@example
P (A + E) P^T = L D L^T
@end example

@end ifinfo
where @math{P} is a permutation matrix, @math{E} is a diagonal
perturbation matrix, @math{L} is unit lower triangular, and
@math{D} is diagonal. If @math{A} is sufficiently positive
definite, then the perturbation matrix @math{E} will be zero
and this method is equivalent to the pivoted Cholesky algorithm.
For indefinite matrices, the perturbation matrix @math{E} is
computed to ensure that @math{A + E} is positive definite and
well conditioned.

@deftypefun int gsl_linalg_mcholesky_decomp (gsl_matrix * @var{A}, gsl_permutation * @var{p}, gsl_vector * @var{E})
This function factors the symmetric, indefinite square matrix
@var{A} into the Modified Cholesky decomposition @math{P (A + E) P^T = L D L^T}.
On input, the values from the diagonal and lower-triangular part of the matrix @var{A} are
used to construct the factorization. On output the diagonal of the input matrix @var{A} stores the diagonal
elements of @math{D}, and the lower triangular portion of @var{A}
contains the matrix @math{L}. Since @math{L} has ones on its diagonal these
do not need to be explicitely stored. The upper triangular portion of @var{A}
is unmodified. The permutation matrix @math{P} is
stored in @var{p} on output. The diagonal perturbation matrix is stored in
@var{E} on output. The parameter @var{E} may be set to NULL if it is not
required.
@end deftypefun

@deftypefun int gsl_linalg_mcholesky_solve (const gsl_matrix * @var{LDLT}, const gsl_permutation * @var{p}, const gsl_vector * @var{b}, gsl_vector * @var{x})
This function solves the perturbed system @math{(A + E) x = b} using the Cholesky
decomposition of @math{A + E} held in the matrix @var{LDLT} and permutation
@var{p} which must have been previously computed by @code{gsl_linalg_mcholesky_decomp}.
@end deftypefun

@deftypefun int gsl_linalg_mcholesky_svx (const gsl_matrix * @var{LDLT}, const gsl_permutation * @var{p}, gsl_vector * @var{x})
This function solves the perturbed system @math{(A + E) x = b} in-place using the Cholesky
decomposition of @math{A + E} held in the matrix @var{LDLT} and permutation
@var{p} which must have been previously computed by @code{gsl_linalg_mcholesky_decomp}.
On input, @var{x} contains the right hand side vector @math{b} which is
replaced by the solution vector on output.
@end deftypefun

@deftypefun int gsl_linalg_mcholesky_rcond (const gsl_matrix * @var{LDLT}, const gsl_permutation * @var{p}, double * @var{rcond}, gsl_vector * @var{work})
This function estimates the reciprocal condition number (using the 1-norm) of the perturbed matrix
@math{A + E}, using its pivoted Cholesky decomposition provided in @var{LDLT}.
The reciprocal condition number estimate, defined as @math{1 / (||A + E||_1 \cdot ||(A + E)^{-1}||_1)}, is stored
in @var{rcond}.  Additional workspace of size @math{3 N} is required in @var{work}.
@end deftypefun

@node Tridiagonal Decomposition of Real Symmetric Matrices
@section Tridiagonal Decomposition of Real Symmetric Matrices
@cindex  tridiagonal decomposition

A symmetric matrix @math{A} can be factorized by similarity
transformations into the form,
@tex
\beforedisplay
$$
A = Q T Q^T
$$
\afterdisplay
@end tex
@ifinfo

@example
A = Q T Q^T
@end example

@end ifinfo
@noindent
where @math{Q} is an orthogonal matrix and @math{T} is a symmetric
tridiagonal matrix.

@deftypefun int gsl_linalg_symmtd_decomp (gsl_matrix * @var{A}, gsl_vector * @var{tau})
This function factorizes the symmetric square matrix @var{A} into the
symmetric tridiagonal decomposition @math{Q T Q^T}.  On output the
diagonal and subdiagonal part of the input matrix @var{A} contain the
tridiagonal matrix @math{T}.  The remaining lower triangular part of the
input matrix contains the Householder vectors which, together with the
Householder coefficients @var{tau}, encode the orthogonal matrix
@math{Q}. This storage scheme is the same as used by @sc{lapack}.  The
upper triangular part of @var{A} is not referenced.
@end deftypefun

@deftypefun int gsl_linalg_symmtd_unpack (const gsl_matrix * @var{A}, const gsl_vector * @var{tau}, gsl_matrix * @var{Q}, gsl_vector * @var{diag}, gsl_vector * @var{subdiag})
This function unpacks the encoded symmetric tridiagonal decomposition
(@var{A}, @var{tau}) obtained from @code{gsl_linalg_symmtd_decomp} into
the orthogonal matrix @var{Q}, the vector of diagonal elements @var{diag}
and the vector of subdiagonal elements @var{subdiag}.  
@end deftypefun

@deftypefun int gsl_linalg_symmtd_unpack_T (const gsl_matrix * @var{A}, gsl_vector * @var{diag}, gsl_vector * @var{subdiag})
This function unpacks the diagonal and subdiagonal of the encoded
symmetric tridiagonal decomposition (@var{A}, @var{tau}) obtained from
@code{gsl_linalg_symmtd_decomp} into the vectors @var{diag} and @var{subdiag}.
@end deftypefun

@node Tridiagonal Decomposition of Hermitian Matrices
@section Tridiagonal Decomposition of Hermitian Matrices
@cindex tridiagonal decomposition

A hermitian matrix @math{A} can be factorized by similarity
transformations into the form,
@tex
\beforedisplay
$$
A = U T U^T
$$
\afterdisplay
@end tex
@ifinfo

@example
A = U T U^T
@end example

@end ifinfo
@noindent
where @math{U} is a unitary matrix and @math{T} is a real symmetric
tridiagonal matrix.


@deftypefun int gsl_linalg_hermtd_decomp (gsl_matrix_complex * @var{A}, gsl_vector_complex * @var{tau})
This function factorizes the hermitian matrix @var{A} into the symmetric
tridiagonal decomposition @math{U T U^T}.  On output the real parts of
the diagonal and subdiagonal part of the input matrix @var{A} contain
the tridiagonal matrix @math{T}.  The remaining lower triangular part of
the input matrix contains the Householder vectors which, together with
the Householder coefficients @var{tau}, encode the unitary matrix
@math{U}. This storage scheme is the same as used by @sc{lapack}.  The
upper triangular part of @var{A} and imaginary parts of the diagonal are
not referenced.
@end deftypefun

@deftypefun int gsl_linalg_hermtd_unpack (const gsl_matrix_complex * @var{A}, const gsl_vector_complex * @var{tau}, gsl_matrix_complex * @var{U}, gsl_vector * @var{diag}, gsl_vector * @var{subdiag})
This function unpacks the encoded tridiagonal decomposition (@var{A},
@var{tau}) obtained from @code{gsl_linalg_hermtd_decomp} into the
unitary matrix @var{U}, the real vector of diagonal elements @var{diag} and
the real vector of subdiagonal elements @var{subdiag}. 
@end deftypefun

@deftypefun int gsl_linalg_hermtd_unpack_T (const gsl_matrix_complex * @var{A}, gsl_vector * @var{diag}, gsl_vector * @var{subdiag})
This function unpacks the diagonal and subdiagonal of the encoded
tridiagonal decomposition (@var{A}, @var{tau}) obtained from the
@code{gsl_linalg_hermtd_decomp} into the real vectors
@var{diag} and @var{subdiag}.
@end deftypefun

@node Hessenberg Decomposition of Real Matrices
@section Hessenberg Decomposition of Real Matrices
@cindex Hessenberg decomposition

A general real matrix @math{A} can be decomposed by orthogonal
similarity transformations into the form
@tex
\beforedisplay
$$
A = U H U^T
$$
\afterdisplay
@end tex
@ifinfo

@example
A = U H U^T
@end example

@end ifinfo
where @math{U} is orthogonal and @math{H} is an upper Hessenberg matrix,
meaning that it has zeros below the first subdiagonal. The
Hessenberg reduction is the first step in the Schur decomposition
for the nonsymmetric eigenvalue problem, but has applications in
other areas as well.

@deftypefun int gsl_linalg_hessenberg_decomp (gsl_matrix * @var{A}, gsl_vector * @var{tau})
This function computes the Hessenberg decomposition of the matrix
@var{A} by applying the similarity transformation @math{H = U^T A U}.
On output, @math{H} is stored in the upper portion of @var{A}. The
information required to construct the matrix @math{U} is stored in
the lower triangular portion of @var{A}. @math{U} is a product
of @math{N - 2} Householder matrices. The Householder vectors
are stored in the lower portion of @var{A} (below the subdiagonal)
and the Householder coefficients are stored in the vector @var{tau}.
@var{tau} must be of length @var{N}.
@end deftypefun

@deftypefun int gsl_linalg_hessenberg_unpack (gsl_matrix * @var{H}, gsl_vector * @var{tau}, gsl_matrix * @var{U})
This function constructs the orthogonal matrix @math{U} from the
information stored in the Hessenberg matrix @var{H} along with the
vector @var{tau}. @var{H} and @var{tau} are outputs from
@code{gsl_linalg_hessenberg_decomp}.
@end deftypefun

@deftypefun int gsl_linalg_hessenberg_unpack_accum (gsl_matrix * @var{H}, gsl_vector * @var{tau}, gsl_matrix * @var{V})
This function is similar to @code{gsl_linalg_hessenberg_unpack}, except
it accumulates the matrix @var{U} into @var{V}, so that @math{V' = VU}.
The matrix @var{V} must be initialized prior to calling this function.
Setting @var{V} to the identity matrix provides the same result as
@code{gsl_linalg_hessenberg_unpack}. If @var{H} is order @var{N}, then
@var{V} must have @var{N} columns but may have any number of rows.
@end deftypefun

@deftypefun int gsl_linalg_hessenberg_set_zero (gsl_matrix * @var{H})
This function sets the lower triangular portion of @var{H}, below
the subdiagonal, to zero. It is useful for clearing out the
Householder vectors after calling @code{gsl_linalg_hessenberg_decomp}.
@end deftypefun

@node Hessenberg-Triangular Decomposition of Real Matrices
@section Hessenberg-Triangular Decomposition of Real Matrices
@cindex Hessenberg triangular decomposition

A general real matrix pair (@math{A}, @math{B}) can be decomposed by
orthogonal similarity transformations into the form
@tex
\beforedisplay
$$
A = U H V^T
$$
$$
B = U R V^T
$$
\afterdisplay
@end tex
@ifinfo

@example
A = U H V^T
B = U R V^T
@end example

@end ifinfo
where @math{U} and @math{V} are orthogonal, @math{H} is an upper
Hessenberg matrix, and @math{R} is upper triangular. The
Hessenberg-Triangular reduction is the first step in the generalized
Schur decomposition for the generalized eigenvalue problem.

@deftypefun int gsl_linalg_hesstri_decomp (gsl_matrix * @var{A}, gsl_matrix * @var{B}, gsl_matrix * @var{U}, gsl_matrix * @var{V}, gsl_vector * @var{work})
This function computes the Hessenberg-Triangular decomposition of the
matrix pair (@var{A}, @var{B}). On output, @math{H} is stored in @var{A},
and @math{R} is stored in @var{B}. If @var{U} and @var{V} are provided
(they may be null), the similarity transformations are stored in them.
Additional workspace of length @math{N} is needed in @var{work}.
@end deftypefun

@node Bidiagonalization
@section Bidiagonalization 
@cindex bidiagonalization of real matrices

A general matrix @math{A} can be factorized by similarity
transformations into the form,
@tex
\beforedisplay
$$
A = U B V^T
$$
\afterdisplay
@end tex
@ifinfo

@example
A = U B V^T
@end example

@end ifinfo
@noindent
where @math{U} and @math{V} are orthogonal matrices and @math{B} is a
@math{N}-by-@math{N} bidiagonal matrix with non-zero entries only on the
diagonal and superdiagonal.  The size of @var{U} is @math{M}-by-@math{N}
and the size of @var{V} is @math{N}-by-@math{N}.

@deftypefun int gsl_linalg_bidiag_decomp (gsl_matrix * @var{A}, gsl_vector * @var{tau_U}, gsl_vector * @var{tau_V})
This function factorizes the @math{M}-by-@math{N} matrix @var{A} into
bidiagonal form @math{U B V^T}.  The diagonal and superdiagonal of the
matrix @math{B} are stored in the diagonal and superdiagonal of @var{A}.
The orthogonal matrices @math{U} and @var{V} are stored as compressed
Householder vectors in the remaining elements of @var{A}.  The
Householder coefficients are stored in the vectors @var{tau_U} and
@var{tau_V}.  The length of @var{tau_U} must equal the number of
elements in the diagonal of @var{A} and the length of @var{tau_V} should
be one element shorter.
@end deftypefun

@deftypefun int gsl_linalg_bidiag_unpack (const gsl_matrix * @var{A}, const gsl_vector * @var{tau_U}, gsl_matrix * @var{U}, const gsl_vector * @var{tau_V}, gsl_matrix * @var{V}, gsl_vector * @var{diag}, gsl_vector * @var{superdiag})
This function unpacks the bidiagonal decomposition of @var{A} produced by
@code{gsl_linalg_bidiag_decomp}, (@var{A}, @var{tau_U}, @var{tau_V})
into the separate orthogonal matrices @var{U}, @var{V} and the diagonal
vector @var{diag} and superdiagonal @var{superdiag}.  Note that @var{U}
is stored as a compact @math{M}-by-@math{N} orthogonal matrix satisfying
@math{U^T U = I} for efficiency.
@end deftypefun

@deftypefun int gsl_linalg_bidiag_unpack2 (gsl_matrix * @var{A}, gsl_vector * @var{tau_U}, gsl_vector * @var{tau_V}, gsl_matrix * @var{V})
This function unpacks the bidiagonal decomposition of @var{A} produced by
@code{gsl_linalg_bidiag_decomp}, (@var{A}, @var{tau_U}, @var{tau_V})
into the separate orthogonal matrices @var{U}, @var{V} and the diagonal
vector @var{diag} and superdiagonal @var{superdiag}.  The matrix @var{U}
is stored in-place in @var{A}.
@end deftypefun

@deftypefun int gsl_linalg_bidiag_unpack_B (const gsl_matrix * @var{A}, gsl_vector * @var{diag}, gsl_vector * @var{superdiag})
This function unpacks the diagonal and superdiagonal of the bidiagonal
decomposition of @var{A} from @code{gsl_linalg_bidiag_decomp}, into
the diagonal vector @var{diag} and superdiagonal vector @var{superdiag}.
@end deftypefun

@node Givens Rotations
@section Givens Rotations
@cindex Givens rotation

A Givens rotation is a rotation in the plane acting on two elements
of a given vector. It can be represented in matrix form as
@tex
\beforedisplay
$$
G(i,j,\theta) =
\left(
\matrix{
1 & \ldots & 0 & \ldots & 0 & \ldots & 0 \cr
\vdots & \ddots & \vdots &  & \vdots & & \vdots \cr
0 & \ldots & \cos{\theta} & \ldots & -\sin{\theta} & \ldots & 0 \cr
\vdots &  & \vdots & \ddots & \vdots & & \vdots \cr
0 & \ldots & \sin{\theta} & \ldots & \cos{\theta} & \ldots & 0 \cr
\vdots &  & \vdots &  & \vdots & \ddots & \vdots \cr
0 & \ldots & 0 & \ldots & 0 & \ldots & 1 \cr
}
\right)
$$
\afterdisplay
@end tex
@ifinfo
@example
@end example
@end ifinfo
where the @math{\cos{\theta}} and @math{\sin{\theta}} appear at
the intersection of the @math{i}th and @math{j}th rows and columns.
When acting on a vector @math{x}, @math{G(i,j,\theta) x} performs
a rotation of the @math{(i,j)} elements of @math{x}. Givens
rotations are typically used to introduce zeros in
vectors, such as during the QR decomposition of a matrix. In this
case, it is typically desired to find @math{c} and @math{s} such that
@tex
\beforedisplay
$$
\left(
\matrix{
c & -s \cr
s & c
}
\right)
\left(
\matrix{
a \cr
b
}
\right) =
\left(
\matrix{
r \cr
0
}
\right)
$$
\afterdisplay
@end tex
@ifinfo
@example
@end example
@end ifinfo
with @math{r = \sqrt{a^2 + b^2}}.

@deftypefun void gsl_linalg_givens (const double @var{a}, const double @var{b}, double * @var{c}, double * @var{s})
This function computes @math{c = \cos{\theta}} and @math{s = \sin{\theta}}
so that the Givens matrix @math{G(\theta)} acting on the
vector @math{(a,b)} produces @math{(r, 0)}, with @math{r = \sqrt{a^2 + b^2}}.
@end deftypefun

@deftypefun void gsl_linalg_givens_gv (gsl_vector * @var{v}, const size_t @var{i}, const size_t @var{j}, const double @var{c}, const double @var{s})
This function applies the Givens rotation defined by
@math{c = \cos{\theta}} and @math{s = \sin{\theta}} to the @var{i}
and @var{j} elements of @var{v}. On output,
@math{(v(i),v(j)) \leftarrow G(\theta) (v(i),v(j))}.
@end deftypefun

@node Householder Transformations
@section Householder Transformations
@cindex Householder matrix
@cindex Householder transformation
@cindex transformation, Householder

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 @math{P} takes the form,
@tex
\beforedisplay
$$
P = I - \tau v v^T
$$
\afterdisplay
@end tex
@ifinfo

@example
P = I - \tau v v^T
@end example

@end ifinfo
@noindent
where @math{v} is a vector (called the @dfn{Householder vector}) and
@math{\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.

@deftypefun double gsl_linalg_householder_transform (gsl_vector * @var{w})
@deftypefunx gsl_complex gsl_linalg_complex_householder_transform (gsl_vector_complex * @var{w})
This function prepares a Householder transformation @math{P = I - \tau v
v^T} which can be used to zero all the elements of the input vector @var{w}
except the first. On output the Householder vector @var{v} is stored in
@var{w} and the scalar @math{\tau} is returned. The householder vector
@var{v} is normalized so that @var{v}[0] = 1, however this 1 is not
stored in the output vector. Instead, @var{w}[0] is set to
the first element of the transformed vector, so that if
@math{u = P w}, @var{w}[0] = @math{u}[0] on output and the remainder
of @math{u} is zero.
@end deftypefun

@deftypefun int gsl_linalg_householder_hm (double @var{tau}, const gsl_vector * @var{v}, gsl_matrix * @var{A})
@deftypefunx int gsl_linalg_complex_householder_hm (gsl_complex @var{tau}, const gsl_vector_complex * @var{v}, gsl_matrix_complex * @var{A})
This function applies the Householder matrix @math{P} defined by the
scalar @var{tau} and the vector @var{v} to the left-hand side of the
matrix @var{A}. On output the result @math{P A} is stored in @var{A}.
@end deftypefun

@deftypefun int gsl_linalg_householder_mh (double @var{tau}, const gsl_vector * @var{v}, gsl_matrix * @var{A})
@deftypefunx int gsl_linalg_complex_householder_mh (gsl_complex @var{tau}, const gsl_vector_complex * @var{v}, gsl_matrix_complex * @var{A})
This function applies the Householder matrix @math{P} defined by the
scalar @var{tau} and the vector @var{v} to the right-hand side of the
matrix @var{A}. On output the result @math{A P} is stored in @var{A}.
@end deftypefun

@deftypefun int gsl_linalg_householder_hv (double @var{tau}, const gsl_vector * @var{v}, gsl_vector * @var{w})
@deftypefunx int gsl_linalg_complex_householder_hv (gsl_complex @var{tau}, const gsl_vector_complex * @var{v}, gsl_vector_complex * @var{w})
This function applies the Householder transformation @math{P} defined by
the scalar @var{tau} and the vector @var{v} to the vector @var{w}.  On
output the result @math{P w} is stored in @var{w}.
@end deftypefun

@comment @deftypefun int gsl_linalg_householder_hm1 (double tau, gsl_matrix * A)
@comment This function applies the Householder transform, defined by the scalar
@comment @var{tau} and the vector @var{v}, to a matrix being build up from the
@comment identity matrix, using the first column of @var{A} as a householder vector.
@comment @end deftypefun

@node Householder solver for linear systems
@section Householder solver for linear systems
@cindex solution of linear system by Householder transformations
@cindex Householder linear solver

@deftypefun int gsl_linalg_HH_solve (gsl_matrix * @var{A}, const gsl_vector * @var{b}, gsl_vector * @var{x})
This function solves the system @math{A x = b} directly using
Householder transformations. On output the solution is stored in @var{x}
and @var{b} is not modified. The matrix @var{A} is destroyed by the
Householder transformations.
@end deftypefun

@deftypefun int gsl_linalg_HH_svx (gsl_matrix * @var{A}, gsl_vector * @var{x})
This function solves the system @math{A x = b} in-place using
Householder transformations.  On input @var{x} should contain the
right-hand side @math{b}, which is replaced by the solution on output.  The
matrix @var{A} is destroyed by the Householder transformations.
@end deftypefun

@node Tridiagonal Systems
@section Tridiagonal Systems
@cindex tridiagonal systems

The functions described in this section efficiently solve symmetric,
non-symmetric and cyclic tridiagonal systems with minimal storage.
Note that the current implementations of these functions use a variant
of Cholesky decomposition, so the tridiagonal matrix must be positive
definite.  For non-positive definite matrices, the functions return
the error code @code{GSL_ESING}.

@deftypefun int gsl_linalg_solve_tridiag (const gsl_vector * @var{diag}, const gsl_vector * @var{e}, const gsl_vector * @var{f}, const gsl_vector * @var{b}, gsl_vector * @var{x})
This function solves the general @math{N}-by-@math{N} system @math{A x =
b} where @var{A} is tridiagonal (@c{$N\geq 2$}
@math{N >= 2}). The super-diagonal and
sub-diagonal vectors @var{e} and @var{f} must be one element shorter
than the diagonal vector @var{diag}.  The form of @var{A} for the 4-by-4
case is shown below,
@tex
\beforedisplay
$$
A = \pmatrix{d_0&e_0&  0& 0\cr
             f_0&d_1&e_1& 0\cr
             0  &f_1&d_2&e_2\cr 
             0  &0  &f_2&d_3\cr}
$$
\afterdisplay
@end tex
@ifinfo

@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 )
@end example
@end ifinfo
@noindent
@end deftypefun

@deftypefun int gsl_linalg_solve_symm_tridiag (const gsl_vector * @var{diag}, const gsl_vector * @var{e}, const gsl_vector * @var{b}, gsl_vector * @var{x})
This function solves the general @math{N}-by-@math{N} system @math{A x =
b} where @var{A} is symmetric tridiagonal (@c{$N\geq 2$}
@math{N >= 2}).  The off-diagonal vector
@var{e} must be one element shorter than the diagonal vector @var{diag}.
The form of @var{A} for the 4-by-4 case is shown below,
@tex
\beforedisplay
$$
A = \pmatrix{d_0&e_0&  0& 0\cr
             e_0&d_1&e_1& 0\cr
             0  &e_1&d_2&e_2\cr 
             0  &0  &e_2&d_3\cr}
$$
\afterdisplay
@end tex
@ifinfo

@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 )
@end example
@end ifinfo
@end deftypefun

@deftypefun int gsl_linalg_solve_cyc_tridiag (const gsl_vector * @var{diag}, const gsl_vector * @var{e}, const gsl_vector * @var{f}, const gsl_vector * @var{b}, gsl_vector * @var{x})
This function solves the general @math{N}-by-@math{N} system @math{A x =
b} where @var{A} is cyclic tridiagonal (@c{$N\geq 3$}
@math{N >= 3}).  The cyclic super-diagonal and
sub-diagonal vectors @var{e} and @var{f} must have the same number of
elements as the diagonal vector @var{diag}.  The form of @var{A} for the
4-by-4 case is shown below,
@tex
\beforedisplay
$$
A = \pmatrix{d_0&e_0& 0 &f_3\cr
             f_0&d_1&e_1& 0 \cr
              0 &f_1&d_2&e_2\cr 
             e_3& 0 &f_2&d_3\cr}
$$
\afterdisplay
@end tex
@ifinfo

@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 )
@end example
@end ifinfo
@end deftypefun


@deftypefun int gsl_linalg_solve_symm_cyc_tridiag (const gsl_vector * @var{diag}, const gsl_vector * @var{e}, const gsl_vector * @var{b}, gsl_vector * @var{x})
This function solves the general @math{N}-by-@math{N} system @math{A x =
b} where @var{A} is symmetric cyclic tridiagonal (@c{$N\geq 3$}
@math{N >= 3}).  The cyclic
off-diagonal vector @var{e} must have the same number of elements as the
diagonal vector @var{diag}.  The form of @var{A} for the 4-by-4 case is
shown below,
@tex
\beforedisplay
$$
A = \pmatrix{d_0&e_0& 0 &e_3\cr
             e_0&d_1&e_1& 0 \cr
              0 &e_1&d_2&e_2\cr 
             e_3& 0 &e_2&d_3\cr}
$$
\afterdisplay
@end tex
@ifinfo

@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 )
@end example
@end ifinfo
@end deftypefun

@node Triangular Systems
@section Triangular Systems
@cindex triangular systems

@deftypefun int gsl_linalg_tri_upper_invert (gsl_matrix * @var{T})
@deftypefunx int gsl_linalg_tri_lower_invert (gsl_matrix * @var{T})
@deftypefunx int gsl_linalg_tri_upper_unit_invert (gsl_matrix * @var{T})
@deftypefunx int gsl_linalg_tri_lower_unit_invert (gsl_matrix * @var{T})
These functions calculate the in-place inverse of the triangular matrix @var{T}. When
the @code{upper} prefix is specified, then the upper triangle of @var{T} is used, and when
the @code{lower} prefix is specified, the lower triangle is used. If the @code{unit}
prefix is specified, then the diagonal elements of the matrix @var{T} are taken as
unity and are not referenced. Otherwise the diagonal elements are used in the inversion.
@end deftypefun

@deftypefun int gsl_linalg_tri_upper_rcond (const gsl_matrix * @var{T}, double * @var{rcond}, gsl_vector * @var{work})
@deftypefunx int gsl_linalg_tri_lower_rcond (const gsl_matrix * @var{T}, double * @var{rcond}, gsl_vector * @var{work})
These functions estimate the reciprocal condition number, in the 1-norm, of the upper or lower
@math{N}-by-@math{N} triangular matrix @var{T}. The reciprocal condition number
is stored in @var{rcond} on output, and is defined by @math{1 / (||T||_1 \cdot ||T^{-1}||_1)}.
Additional workspace of size @math{3 N} is required in @var{work}.
@end deftypefun

@node Balancing
@section Balancing
@cindex balancing matrices

The process of balancing a matrix applies similarity transformations
to make the rows and columns have comparable norms. This is
useful, for example, to reduce roundoff errors in the solution
of eigenvalue problems. Balancing a matrix @math{A} consists
of replacing @math{A} with a similar matrix
@tex
\beforedisplay
$$
A' = D^{-1} A D
$$
\afterdisplay
@end tex
@ifinfo

@example
A' = D^(-1) A D
@end example

@end ifinfo
where @math{D} is a diagonal matrix whose entries are powers
of the floating point radix.

@deftypefun int gsl_linalg_balance_matrix (gsl_matrix * @var{A}, gsl_vector * @var{D})
This function replaces the matrix @var{A} with its balanced counterpart
and stores the diagonal elements of the similarity transformation
into the vector @var{D}.
@end deftypefun

@node Linear Algebra Examples
@section Examples

The following program solves the linear system @math{A x = b}. The
system to be solved is,
@tex
\beforedisplay
$$
\left(
\matrix{0.18& 0.60& 0.57& 0.96\cr
0.41& 0.24& 0.99& 0.58\cr
0.14& 0.30& 0.97& 0.66\cr
0.51& 0.13& 0.19& 0.85}
\right)
\left(
\matrix{x_0\cr
x_1\cr
x_2\cr
x_3}
\right)
=
\left(
\matrix{1.0\cr
2.0\cr
3.0\cr
4.0}
\right)
$$
\afterdisplay
@end tex
@ifinfo

@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]
@end example

@end ifinfo
@noindent
and the solution is found using LU decomposition of the matrix @math{A}.

@example
@verbatiminclude examples/linalglu.c
@end example

@noindent
Here is the output from the program,

@example
@verbatiminclude examples/linalglu.txt
@end example

@noindent
This can be verified by multiplying the solution @math{x} by the
original matrix @math{A} using @sc{gnu octave},

@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
@end example

@noindent
This reproduces the original right-hand side vector, @math{b}, in
accordance with the equation @math{A x = b}.

@node Linear Algebra References and Further Reading
@section References and Further Reading

Further information on the algorithms described in this section can be
found in the following book,

@itemize @w{}
@item
G. H. Golub, C. F. Van Loan, @cite{Matrix Computations} (3rd Ed, 1996),
Johns Hopkins University Press, ISBN 0-8018-5414-8.
@end itemize

@noindent
The @sc{lapack} library is described in the following manual,

@itemize @w{}
@item
@cite{LAPACK Users' Guide} (Third Edition, 1999), Published by SIAM,
ISBN 0-89871-447-8.

@uref{http://www.netlib.org/lapack} 
@end itemize

@noindent
The @sc{lapack} source code can be found at the website above, along
with an online copy of the users guide.

@noindent
The Modified Golub-Reinsch algorithm is described in the following paper,

@itemize @w{}
@item
T.F. Chan, ``An Improved Algorithm for Computing the Singular Value
Decomposition'', @cite{ACM Transactions on Mathematical Software}, 8
(1982), pp 72--83.
@end itemize

@noindent
The Jacobi algorithm for singular value decomposition is described in
the following papers,

@itemize @w{}
@item
J.C. Nash, ``A one-sided transformation method for the singular value
decomposition and algebraic eigenproblem'', @cite{Computer Journal},
Volume 18, Number 1 (1975), p 74--76

@item
J.C. Nash and S. Shlien ``Simple algorithms for the partial singular
value decomposition'', @cite{Computer Journal}, Volume 30 (1987), p
268--275.

@item
James Demmel, Kre@v{s}imir Veseli@'c, ``Jacobi's Method is more accurate than
QR'', @cite{Lapack Working Note 15} (LAWN-15), October 1989. Available
from netlib, @uref{http://www.netlib.org/lapack/} in the @code{lawns} or
@code{lawnspdf} directories.
@end itemize

@noindent
The algorithm for estimating a matrix condition number is described in
the following paper,

@itemize @w{}
@item
N. J. Higham, "FORTRAN codes for estimating the one-norm of
a real or complex matrix, with applications to condition estimation",
ACM Trans. Math. Soft., vol. 14, no. 4, pp. 381-396, December 1988.
@end itemize