Added script front-end for primer-design code
[htsworkflow.git] / htswanalysis / MACS / lib / gsl / gsl-1.11 / doc / linalg.texi
1 @cindex linear algebra
2 @cindex solution of linear systems, Ax=b
3 @cindex matrix factorization
4 @cindex factorization of matrices
5
6 This chapter describes functions for solving linear systems.  The
7 library provides linear algebra operations which operate directly on
8 the @code{gsl_vector} and @code{gsl_matrix} objects.  These routines
9 use the standard algorithms from Golub & Van Loan's @cite{Matrix
10 Computations} with Level-1 and Level-2 BLAS calls for efficiency.
11
12 The functions described in this chapter are declared in the header file
13 @file{gsl_linalg.h}.
14
15
16 @menu
17 * LU Decomposition::            
18 * QR Decomposition::            
19 * QR Decomposition with Column Pivoting::  
20 * Singular Value Decomposition::  
21 * Cholesky Decomposition::      
22 * Tridiagonal Decomposition of Real Symmetric Matrices::  
23 * Tridiagonal Decomposition of Hermitian Matrices::  
24 * Hessenberg Decomposition of Real Matrices::
25 * Hessenberg-Triangular Decomposition of Real Matrices::
26 * Bidiagonalization::           
27 * Householder Transformations::  
28 * Householder solver for linear systems::  
29 * Tridiagonal Systems::         
30 * Balancing::
31 * Linear Algebra Examples::     
32 * Linear Algebra References and Further Reading::  
33 @end menu
34
35 @node LU Decomposition
36 @section LU Decomposition
37 @cindex LU decomposition
38
39 A general square matrix @math{A} has an @math{LU} decomposition into
40 upper and lower triangular matrices,
41 @tex
42 \beforedisplay
43 $$
44 P A = L U
45 $$
46 \afterdisplay
47 @end tex
48 @ifinfo
49
50 @example
51 P A = L U
52 @end example
53
54 @end ifinfo
55 @noindent
56 where @math{P} is a permutation matrix, @math{L} is unit lower
57 triangular matrix and @math{U} is upper triangular matrix. For square
58 matrices this decomposition can be used to convert the linear system
59 @math{A x = b} into a pair of triangular systems (@math{L y = P b},
60 @math{U x = y}), which can be solved by forward and back-substitution.
61 Note that the @math{LU} decomposition is valid for singular matrices.
62
63 @deftypefun int gsl_linalg_LU_decomp (gsl_matrix * @var{A}, gsl_permutation * @var{p}, int * @var{signum})
64 @deftypefunx int gsl_linalg_complex_LU_decomp (gsl_matrix_complex * @var{A}, gsl_permutation * @var{p}, int * @var{signum})
65 These functions factorize the square matrix @var{A} into the @math{LU}
66 decomposition @math{PA = LU}.  On output the diagonal and upper
67 triangular part of the input matrix @var{A} contain the matrix
68 @math{U}. The lower triangular part of the input matrix (excluding the
69 diagonal) contains @math{L}.  The diagonal elements of @math{L} are
70 unity, and are not stored.
71
72 The permutation matrix @math{P} is encoded in the permutation
73 @var{p}. The @math{j}-th column of the matrix @math{P} is given by the
74 @math{k}-th column of the identity matrix, where @math{k = p_j} the
75 @math{j}-th element of the permutation vector. The sign of the
76 permutation is given by @var{signum}. It has the value @math{(-1)^n},
77 where @math{n} is the number of interchanges in the permutation.
78
79 The algorithm used in the decomposition is Gaussian Elimination with
80 partial pivoting (Golub & Van Loan, @cite{Matrix Computations},
81 Algorithm 3.4.1).
82 @end deftypefun
83
84 @cindex linear systems, solution of
85 @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})
86 @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})
87 These functions solve the square system @math{A x = b} using the @math{LU}
88 decomposition of @math{A} into (@var{LU}, @var{p}) given by
89 @code{gsl_linalg_LU_decomp} or @code{gsl_linalg_complex_LU_decomp}.
90 @end deftypefun
91
92 @deftypefun int gsl_linalg_LU_svx (const gsl_matrix * @var{LU}, const gsl_permutation * @var{p}, gsl_vector * @var{x})
93 @deftypefunx int gsl_linalg_complex_LU_svx (const gsl_matrix_complex * @var{LU}, const gsl_permutation * @var{p}, gsl_vector_complex * @var{x})
94 These functions solve the square system @math{A x = b} in-place using the
95 @math{LU} decomposition of @math{A} into (@var{LU},@var{p}). On input
96 @var{x} should contain the right-hand side @math{b}, which is replaced
97 by the solution on output.
98 @end deftypefun
99
100 @cindex refinement of solutions in linear systems
101 @cindex iterative refinement of solutions in linear systems
102 @cindex linear systems, refinement of solutions
103 @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{residual})
104 @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{residual})
105 These functions apply an iterative improvement to @var{x}, the solution
106 of @math{A x = b}, using the @math{LU} decomposition of @math{A} into
107 (@var{LU},@var{p}). The initial residual @math{r = A x - b} is also
108 computed and stored in @var{residual}.
109 @end deftypefun
110
111 @cindex inverse of a matrix, by LU decomposition
112 @cindex matrix inverse
113 @deftypefun int gsl_linalg_LU_invert (const gsl_matrix * @var{LU}, const gsl_permutation * @var{p}, gsl_matrix * @var{inverse})
114 @deftypefunx int gsl_linalg_complex_LU_invert (const gsl_matrix_complex * @var{LU}, const gsl_permutation * @var{p}, gsl_matrix_complex * @var{inverse})
115 These functions compute the inverse of a matrix @math{A} from its
116 @math{LU} decomposition (@var{LU},@var{p}), storing the result in the
117 matrix @var{inverse}. The inverse is computed by solving the system
118 @math{A x = b} for each column of the identity matrix.  It is preferable
119 to avoid direct use of the inverse whenever possible, as the linear
120 solver functions can obtain the same result more efficiently and
121 reliably (consult any introductory textbook on numerical linear algebra
122 for details).
123 @end deftypefun
124
125 @cindex determinant of a matrix, by LU decomposition
126 @cindex matrix determinant
127 @deftypefun double gsl_linalg_LU_det (gsl_matrix * @var{LU}, int @var{signum})
128 @deftypefunx gsl_complex gsl_linalg_complex_LU_det (gsl_matrix_complex * @var{LU}, int @var{signum})
129 These functions compute the determinant of a matrix @math{A} from its
130 @math{LU} decomposition, @var{LU}. The determinant is computed as the
131 product of the diagonal elements of @math{U} and the sign of the row
132 permutation @var{signum}.
133 @end deftypefun
134
135 @cindex logarithm of the determinant of a matrix
136 @deftypefun double gsl_linalg_LU_lndet (gsl_matrix * @var{LU})
137 @deftypefunx double gsl_linalg_complex_LU_lndet (gsl_matrix_complex * @var{LU})
138 These functions compute the logarithm of the absolute value of the
139 determinant of a matrix @math{A}, @math{\ln|\det(A)|}, from its @math{LU}
140 decomposition, @var{LU}.  This function may be useful if the direct
141 computation of the determinant would overflow or underflow.
142 @end deftypefun
143
144 @cindex sign of the determinant of a matrix
145 @deftypefun int gsl_linalg_LU_sgndet (gsl_matrix * @var{LU}, int @var{signum})
146 @deftypefunx gsl_complex gsl_linalg_complex_LU_sgndet (gsl_matrix_complex * @var{LU}, int @var{signum})
147 These functions compute the sign or phase factor of the determinant of a
148 matrix @math{A}, @math{\det(A)/|\det(A)|}, from its @math{LU} decomposition,
149 @var{LU}.
150 @end deftypefun
151
152 @node QR Decomposition
153 @section QR Decomposition
154 @cindex QR decomposition
155
156 A general rectangular @math{M}-by-@math{N} matrix @math{A} has a
157 @math{QR} decomposition into the product of an orthogonal
158 @math{M}-by-@math{M} square matrix @math{Q} (where @math{Q^T Q = I}) and
159 an @math{M}-by-@math{N} right-triangular matrix @math{R},
160 @tex
161 \beforedisplay
162 $$
163 A = Q R
164 $$
165 \afterdisplay
166 @end tex
167 @ifinfo
168
169 @example
170 A = Q R
171 @end example
172
173 @end ifinfo
174 @noindent
175 This decomposition can be used to convert the linear system @math{A x =
176 b} into the triangular system @math{R x = Q^T b}, which can be solved by
177 back-substitution. Another use of the @math{QR} decomposition is to
178 compute an orthonormal basis for a set of vectors. The first @math{N}
179 columns of @math{Q} form an orthonormal basis for the range of @math{A},
180 @math{ran(A)}, when @math{A} has full column rank.
181
182 @deftypefun int gsl_linalg_QR_decomp (gsl_matrix * @var{A}, gsl_vector * @var{tau})
183 This function factorizes the @math{M}-by-@math{N} matrix @var{A} into
184 the @math{QR} decomposition @math{A = Q R}.  On output the diagonal and
185 upper triangular part of the input matrix contain the matrix
186 @math{R}. The vector @var{tau} and the columns of the lower triangular
187 part of the matrix @var{A} contain the Householder coefficients and
188 Householder vectors which encode the orthogonal matrix @var{Q}.  The
189 vector @var{tau} must be of length @math{k=\min(M,N)}. The matrix
190 @math{Q} is related to these components by, @math{Q = Q_k ... Q_2 Q_1}
191 where @math{Q_i = I - \tau_i v_i v_i^T} and @math{v_i} is the
192 Householder vector @math{v_i =
193 (0,...,1,A(i+1,i),A(i+2,i),...,A(m,i))}. This is the same storage scheme
194 as used by @sc{lapack}.
195
196 The algorithm used to perform the decomposition is Householder QR (Golub
197 & Van Loan, @cite{Matrix Computations}, Algorithm 5.2.1).
198 @end deftypefun
199
200 @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})
201 This function solves the square system @math{A x = b} using the @math{QR}
202 decomposition of @math{A} into (@var{QR}, @var{tau}) given by
203 @code{gsl_linalg_QR_decomp}. The least-squares solution for rectangular systems can
204 be found using @code{gsl_linalg_QR_lssolve}.
205 @end deftypefun
206
207 @deftypefun int gsl_linalg_QR_svx (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, gsl_vector * @var{x})
208 This function solves the square system @math{A x = b} in-place using the
209 @math{QR} decomposition of @math{A} into (@var{QR},@var{tau}) given by
210 @code{gsl_linalg_QR_decomp}. On input @var{x} should contain the
211 right-hand side @math{b}, which is replaced by the solution on output.
212 @end deftypefun
213
214 @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})
215 This function finds the least squares solution to the overdetermined
216 system @math{A x = b} where the matrix @var{A} has more rows than
217 columns.  The least squares solution minimizes the Euclidean norm of the
218 residual, @math{||Ax - b||}.The routine uses the @math{QR} decomposition
219 of @math{A} into (@var{QR}, @var{tau}) given by
220 @code{gsl_linalg_QR_decomp}.  The solution is returned in @var{x}.  The
221 residual is computed as a by-product and stored in @var{residual}.
222 @end deftypefun
223
224 @deftypefun int gsl_linalg_QR_QTvec (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, gsl_vector * @var{v})
225 This function applies the matrix @math{Q^T} encoded in the decomposition
226 (@var{QR},@var{tau}) to the vector @var{v}, storing the result @math{Q^T
227 v} in @var{v}.  The matrix multiplication is carried out directly using
228 the encoding of the Householder vectors without needing to form the full
229 matrix @math{Q^T}.
230 @end deftypefun
231
232 @deftypefun int gsl_linalg_QR_Qvec (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, gsl_vector * @var{v})
233 This function applies the matrix @math{Q} encoded in the decomposition
234 (@var{QR},@var{tau}) to the vector @var{v}, storing the result @math{Q
235 v} in @var{v}.  The matrix multiplication is carried out directly using
236 the encoding of the Householder vectors without needing to form the full
237 matrix @math{Q}.
238 @end deftypefun
239
240 @deftypefun int gsl_linalg_QR_QTmat (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, gsl_matrix * @var{A})
241 This function applies the matrix @math{Q^T} encoded in the decomposition
242 (@var{QR},@var{tau}) to the matrix @var{A}, storing the result @math{Q^T
243 A} in @var{A}.  The matrix multiplication is carried out directly using
244 the encoding of the Householder vectors without needing to form the full
245 matrix @math{Q^T}.
246 @end deftypefun
247
248 @deftypefun int gsl_linalg_QR_Rsolve (const gsl_matrix * @var{QR}, const gsl_vector * @var{b}, gsl_vector * @var{x})
249 This function solves the triangular system @math{R x = b} for
250 @var{x}. It may be useful if the product @math{b' = Q^T b} has already
251 been computed using @code{gsl_linalg_QR_QTvec}.
252 @end deftypefun
253
254 @deftypefun int gsl_linalg_QR_Rsvx (const gsl_matrix * @var{QR}, gsl_vector * @var{x})
255 This function solves the triangular system @math{R x = b} for @var{x}
256 in-place. On input @var{x} should contain the right-hand side @math{b}
257 and is replaced by the solution on output. This function may be useful if
258 the product @math{b' = Q^T b} has already been computed using
259 @code{gsl_linalg_QR_QTvec}.
260 @end deftypefun
261
262 @deftypefun int gsl_linalg_QR_unpack (const gsl_matrix * @var{QR}, const gsl_vector * @var{tau}, gsl_matrix * @var{Q}, gsl_matrix * @var{R})
263 This function unpacks the encoded @math{QR} decomposition
264 (@var{QR},@var{tau}) into the matrices @var{Q} and @var{R}, where
265 @var{Q} is @math{M}-by-@math{M} and @var{R} is @math{M}-by-@math{N}. 
266 @end deftypefun
267
268 @deftypefun int gsl_linalg_QR_QRsolve (gsl_matrix * @var{Q}, gsl_matrix * @var{R}, const gsl_vector * @var{b}, gsl_vector * @var{x})
269 This function solves the system @math{R x = Q^T b} for @var{x}. It can
270 be used when the @math{QR} decomposition of a matrix is available in
271 unpacked form as (@var{Q}, @var{R}).
272 @end deftypefun
273
274 @deftypefun int gsl_linalg_QR_update (gsl_matrix * @var{Q}, gsl_matrix * @var{R}, gsl_vector * @var{w}, const gsl_vector * @var{v})
275 This function performs a rank-1 update @math{w v^T} of the @math{QR}
276 decomposition (@var{Q}, @var{R}). The update is given by @math{Q'R' = Q
277 R + w v^T} where the output matrices @math{Q'} and @math{R'} are also
278 orthogonal and right triangular. Note that @var{w} is destroyed by the
279 update.
280 @end deftypefun
281
282 @deftypefun int gsl_linalg_R_solve (const gsl_matrix * @var{R}, const gsl_vector * @var{b}, gsl_vector * @var{x})
283 This function solves the triangular system @math{R x = b} for the
284 @math{N}-by-@math{N} matrix @var{R}.
285 @end deftypefun
286
287 @deftypefun int gsl_linalg_R_svx (const gsl_matrix * @var{R}, gsl_vector * @var{x})
288 This function solves the triangular system @math{R x = b} in-place. On
289 input @var{x} should contain the right-hand side @math{b}, which is
290 replaced by the solution on output.
291 @end deftypefun
292
293 @node QR Decomposition with Column Pivoting
294 @section QR Decomposition with Column Pivoting
295 @cindex QR decomposition with column pivoting
296
297 The @math{QR} decomposition can be extended to the rank deficient case
298 by introducing a column permutation @math{P},
299 @tex
300 \beforedisplay
301 $$
302 A P = Q R
303 $$
304 \afterdisplay
305 @end tex
306 @ifinfo
307
308 @example
309 A P = Q R
310 @end example
311
312 @end ifinfo
313 @noindent
314 The first @math{r} columns of @math{Q} form an orthonormal basis
315 for the range of @math{A} for a matrix with column rank @math{r}.  This
316 decomposition can also be used to convert the linear system @math{A x =
317 b} into the triangular system @math{R y = Q^T b, x = P y}, which can be
318 solved by back-substitution and permutation.  We denote the @math{QR}
319 decomposition with column pivoting by @math{QRP^T} since @math{A = Q R
320 P^T}.
321
322 @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})
323 This function factorizes the @math{M}-by-@math{N} matrix @var{A} into
324 the @math{QRP^T} decomposition @math{A = Q R P^T}.  On output the
325 diagonal and upper triangular part of the input matrix contain the
326 matrix @math{R}. The permutation matrix @math{P} is stored in the
327 permutation @var{p}.  The sign of the permutation is given by
328 @var{signum}. It has the value @math{(-1)^n}, where @math{n} is the
329 number of interchanges in the permutation. The vector @var{tau} and the
330 columns of the lower triangular part of the matrix @var{A} contain the
331 Householder coefficients and vectors which encode the orthogonal matrix
332 @var{Q}.  The vector @var{tau} must be of length @math{k=\min(M,N)}. The
333 matrix @math{Q} is related to these components by, @math{Q = Q_k ... Q_2
334 Q_1} where @math{Q_i = I - \tau_i v_i v_i^T} and @math{v_i} is the
335 Householder vector @math{v_i =
336 (0,...,1,A(i+1,i),A(i+2,i),...,A(m,i))}. This is the same storage scheme
337 as used by @sc{lapack}.  The vector @var{norm} is a workspace of length
338 @var{N} used for column pivoting.
339
340 The algorithm used to perform the decomposition is Householder QR with
341 column pivoting (Golub & Van Loan, @cite{Matrix Computations}, Algorithm
342 5.4.1).
343 @end deftypefun
344
345 @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})
346 This function factorizes the matrix @var{A} into the decomposition
347 @math{A = Q R P^T} without modifying @var{A} itself and storing the
348 output in the separate matrices @var{q} and @var{r}.
349 @end deftypefun
350
351 @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})
352 This function solves the square system @math{A x = b} using the @math{QRP^T}
353 decomposition of @math{A} into (@var{QR}, @var{tau}, @var{p}) given by
354 @code{gsl_linalg_QRPT_decomp}.
355 @end deftypefun
356
357 @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})
358 This function solves the square system @math{A x = b} in-place using the
359 @math{QRP^T} decomposition of @math{A} into
360 (@var{QR},@var{tau},@var{p}). On input @var{x} should contain the
361 right-hand side @math{b}, which is replaced by the solution on output.
362 @end deftypefun
363
364 @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})
365 This function solves the square system @math{R P^T x = Q^T b} for
366 @var{x}. It can be used when the @math{QR} decomposition of a matrix is
367 available in unpacked form as (@var{Q}, @var{R}).
368 @end deftypefun
369
370 @deftypefun int gsl_linalg_QRPT_update (gsl_matrix * @var{Q}, gsl_matrix * @var{R}, const gsl_permutation * @var{p}, gsl_vector * @var{u}, const gsl_vector * @var{v})
371 This function performs a rank-1 update @math{w v^T} of the @math{QRP^T}
372 decomposition (@var{Q}, @var{R}, @var{p}). The update is given by
373 @math{Q'R' = Q R + w v^T} where the output matrices @math{Q'} and
374 @math{R'} are also orthogonal and right triangular. Note that @var{w} is
375 destroyed by the update. The permutation @var{p} is not changed.
376 @end deftypefun
377
378 @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})
379 This function solves the triangular system @math{R P^T x = b} for the
380 @math{N}-by-@math{N} matrix @math{R} contained in @var{QR}.
381 @end deftypefun
382
383 @deftypefun int gsl_linalg_QRPT_Rsvx (const gsl_matrix * @var{QR}, const gsl_permutation * @var{p}, gsl_vector * @var{x})
384 This function solves the triangular system @math{R P^T x = b} in-place
385 for the @math{N}-by-@math{N} matrix @math{R} contained in @var{QR}. On
386 input @var{x} should contain the right-hand side @math{b}, which is
387 replaced by the solution on output.
388 @end deftypefun
389
390 @node Singular Value Decomposition
391 @section Singular Value Decomposition
392 @cindex SVD
393 @cindex singular value decomposition
394
395 A general rectangular @math{M}-by-@math{N} matrix @math{A} has a
396 singular value decomposition (@sc{svd}) into the product of an
397 @math{M}-by-@math{N} orthogonal matrix @math{U}, an @math{N}-by-@math{N}
398 diagonal matrix of singular values @math{S} and the transpose of an
399 @math{N}-by-@math{N} orthogonal square matrix @math{V},
400 @tex
401 \beforedisplay
402 $$
403 A = U S V^T
404 $$
405 \afterdisplay
406 @end tex
407 @ifinfo
408
409 @example
410 A = U S V^T
411 @end example
412
413 @end ifinfo
414 @noindent
415 The singular values
416 @c{$\sigma_i = S_{ii}$}
417 @math{\sigma_i = S_@{ii@}} are all non-negative and are
418 generally chosen to form a non-increasing sequence 
419 @c{$\sigma_1 \ge \sigma_2 \ge ... \ge \sigma_N \ge 0$}
420 @math{\sigma_1 >= \sigma_2 >= ... >= \sigma_N >= 0}.
421
422 The singular value decomposition of a matrix has many practical uses.
423 The condition number of the matrix is given by the ratio of the largest
424 singular value to the smallest singular value. The presence of a zero
425 singular value indicates that the matrix is singular. The number of
426 non-zero singular values indicates the rank of the matrix.  In practice
427 singular value decomposition of a rank-deficient matrix will not produce
428 exact zeroes for singular values, due to finite numerical
429 precision.  Small singular values should be edited by choosing a suitable
430 tolerance.
431
432 For a rank-deficient matrix, the null space of @math{A} is given by
433 the columns of @math{V} corresponding to the zero singular values.
434 Similarly, the range of @math{A} is given by columns of @math{U}
435 corresponding to the non-zero singular values.
436
437 @deftypefun int gsl_linalg_SV_decomp (gsl_matrix * @var{A}, gsl_matrix * @var{V}, gsl_vector * @var{S}, gsl_vector * @var{work})
438 This function factorizes the @math{M}-by-@math{N} matrix @var{A} into
439 the singular value decomposition @math{A = U S V^T} for @c{$M \ge N$}
440 @math{M >= N}.  On output the matrix @var{A} is replaced by
441 @math{U}. The diagonal elements of the singular value matrix @math{S}
442 are stored in the vector @var{S}. The singular values are non-negative
443 and form a non-increasing sequence from @math{S_1} to @math{S_N}. The
444 matrix @var{V} contains the elements of @math{V} in untransposed
445 form. To form the product @math{U S V^T} it is necessary to take the
446 transpose of @var{V}.  A workspace of length @var{N} is required in
447 @var{work}.
448
449 This routine uses the Golub-Reinsch SVD algorithm.
450 @end deftypefun
451
452 @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})
453 This function computes the SVD using the modified Golub-Reinsch
454 algorithm, which is faster for @c{$M \gg N$}
455 @math{M>>N}.  It requires the vector @var{work} of length @var{N} and the
456 @math{N}-by-@math{N} matrix @var{X} as additional working space.
457 @end deftypefun
458
459 @deftypefun int gsl_linalg_SV_decomp_jacobi (gsl_matrix * @var{A}, gsl_matrix * @var{V}, gsl_vector * @var{S})
460 This function computes the SVD of the @math{M}-by-@math{N} matrix @var{A}
461 using one-sided Jacobi orthogonalization for @c{$M \ge N$} 
462 @math{M >= N}.  The Jacobi method can compute singular values to higher
463 relative accuracy than Golub-Reinsch algorithms (see references for
464 details).
465 @end deftypefun
466
467 @deftypefun int gsl_linalg_SV_solve (gsl_matrix * @var{U}, gsl_matrix * @var{V}, gsl_vector * @var{S}, const gsl_vector * @var{b}, gsl_vector * @var{x})
468 This function solves the system @math{A x = b} using the singular value
469 decomposition (@var{U}, @var{S}, @var{V}) of @math{A} given by
470 @code{gsl_linalg_SV_decomp}.
471
472 Only non-zero singular values are used in computing the solution. The
473 parts of the solution corresponding to singular values of zero are
474 ignored.  Other singular values can be edited out by setting them to
475 zero before calling this function. 
476
477 In the over-determined case where @var{A} has more rows than columns the
478 system is solved in the least squares sense, returning the solution
479 @var{x} which minimizes @math{||A x - b||_2}.
480 @end deftypefun
481
482 @node Cholesky Decomposition
483 @section Cholesky Decomposition
484 @cindex Cholesky decomposition
485 @cindex square root of a matrix, Cholesky decomposition
486 @cindex matrix square root, Cholesky decomposition
487
488 A symmetric, positive definite square matrix @math{A} has a Cholesky
489 decomposition into a product of a lower triangular matrix @math{L} and
490 its transpose @math{L^T},
491 @tex
492 \beforedisplay
493 $$
494 A = L L^T
495 $$
496 \afterdisplay
497 @end tex
498 @ifinfo
499
500 @example
501 A = L L^T
502 @end example
503
504 @end ifinfo
505 @noindent
506 This is sometimes referred to as taking the square-root of a matrix. The
507 Cholesky decomposition can only be carried out when all the eigenvalues
508 of the matrix are positive.  This decomposition can be used to convert
509 the linear system @math{A x = b} into a pair of triangular systems
510 (@math{L y = b}, @math{L^T x = y}), which can be solved by forward and
511 back-substitution.
512
513 @deftypefun int gsl_linalg_cholesky_decomp (gsl_matrix * @var{A})
514 @deftypefunx int gsl_linalg_complex_cholesky_decomp (gsl_matrix_complex * @var{A})
515 These functions factorize the symmetric, positive-definite square matrix
516 @var{A} into the Cholesky decomposition @math{A = L L^T} (or
517 @c{$A = L L^{\dagger}$}
518 @math{A = L L^H}
519 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
520 matrix @var{A} contain the matrix @math{L}, while the upper triangular part
521 of the input matrix is overwritten with @math{L^T} (the diagonal terms being
522 identical for both @math{L} and @math{L^T}).  If the matrix is not
523 positive-definite then the decomposition will fail, returning the
524 error code @code{GSL_EDOM}.  
525
526 When testing whether a matrix is positive-definite, disable the error
527 handler first to avoid triggering an error.
528 @end deftypefun
529
530 @deftypefun int gsl_linalg_cholesky_solve (const gsl_matrix * @var{cholesky}, const gsl_vector * @var{b}, gsl_vector * @var{x})
531 @deftypefunx int gsl_linalg_complex_cholesky_solve (const gsl_matrix_complex * @var{cholesky}, const gsl_vector_complex * @var{b}, gsl_vector_complex * @var{x})
532 These functions solve the system @math{A x = b} using the Cholesky
533 decomposition of @math{A} into the matrix @var{cholesky} given by
534 @code{gsl_linalg_cholesky_decomp} or
535 @code{gsl_linalg_complex_cholesky_decomp}.
536 @end deftypefun
537
538 @deftypefun int gsl_linalg_cholesky_svx (const gsl_matrix * @var{cholesky}, gsl_vector * @var{x})
539 @deftypefunx int gsl_linalg_complex_cholesky_svx (const gsl_matrix_complex * @var{cholesky}, gsl_vector_complex * @var{x})
540 These functions solve the system @math{A x = b} in-place using the
541 Cholesky decomposition of @math{A} into the matrix @var{cholesky} given
542 by @code{gsl_linalg_cholesky_decomp} or
543 @code{gsl_linalg_complex_cholesky_decomp}. On input @var{x} should contain
544 the right-hand side @math{b}, which is replaced by the solution on
545 output.
546 @end deftypefun
547
548 @node Tridiagonal Decomposition of Real Symmetric Matrices
549 @section Tridiagonal Decomposition of Real Symmetric Matrices
550 @cindex  tridiagonal decomposition
551
552 A symmetric matrix @math{A} can be factorized by similarity
553 transformations into the form,
554 @tex
555 \beforedisplay
556 $$
557 A = Q T Q^T
558 $$
559 \afterdisplay
560 @end tex
561 @ifinfo
562
563 @example
564 A = Q T Q^T
565 @end example
566
567 @end ifinfo
568 @noindent
569 where @math{Q} is an orthogonal matrix and @math{T} is a symmetric
570 tridiagonal matrix.
571
572 @deftypefun int gsl_linalg_symmtd_decomp (gsl_matrix * @var{A}, gsl_vector * @var{tau})
573 This function factorizes the symmetric square matrix @var{A} into the
574 symmetric tridiagonal decomposition @math{Q T Q^T}.  On output the
575 diagonal and subdiagonal part of the input matrix @var{A} contain the
576 tridiagonal matrix @math{T}.  The remaining lower triangular part of the
577 input matrix contains the Householder vectors which, together with the
578 Householder coefficients @var{tau}, encode the orthogonal matrix
579 @math{Q}. This storage scheme is the same as used by @sc{lapack}.  The
580 upper triangular part of @var{A} is not referenced.
581 @end deftypefun
582
583 @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})
584 This function unpacks the encoded symmetric tridiagonal decomposition
585 (@var{A}, @var{tau}) obtained from @code{gsl_linalg_symmtd_decomp} into
586 the orthogonal matrix @var{Q}, the vector of diagonal elements @var{diag}
587 and the vector of subdiagonal elements @var{subdiag}.  
588 @end deftypefun
589
590 @deftypefun int gsl_linalg_symmtd_unpack_T (const gsl_matrix * @var{A}, gsl_vector * @var{diag}, gsl_vector * @var{subdiag})
591 This function unpacks the diagonal and subdiagonal of the encoded
592 symmetric tridiagonal decomposition (@var{A}, @var{tau}) obtained from
593 @code{gsl_linalg_symmtd_decomp} into the vectors @var{diag} and @var{subdiag}.
594 @end deftypefun
595
596 @node Tridiagonal Decomposition of Hermitian Matrices
597 @section Tridiagonal Decomposition of Hermitian Matrices
598 @cindex tridiagonal decomposition
599
600 A hermitian matrix @math{A} can be factorized by similarity
601 transformations into the form,
602 @tex
603 \beforedisplay
604 $$
605 A = U T U^T
606 $$
607 \afterdisplay
608 @end tex
609 @ifinfo
610
611 @example
612 A = U T U^T
613 @end example
614
615 @end ifinfo
616 @noindent
617 where @math{U} is a unitary matrix and @math{T} is a real symmetric
618 tridiagonal matrix.
619
620
621 @deftypefun int gsl_linalg_hermtd_decomp (gsl_matrix_complex * @var{A}, gsl_vector_complex * @var{tau})
622 This function factorizes the hermitian matrix @var{A} into the symmetric
623 tridiagonal decomposition @math{U T U^T}.  On output the real parts of
624 the diagonal and subdiagonal part of the input matrix @var{A} contain
625 the tridiagonal matrix @math{T}.  The remaining lower triangular part of
626 the input matrix contains the Householder vectors which, together with
627 the Householder coefficients @var{tau}, encode the orthogonal matrix
628 @math{Q}. This storage scheme is the same as used by @sc{lapack}.  The
629 upper triangular part of @var{A} and imaginary parts of the diagonal are
630 not referenced.
631 @end deftypefun
632
633 @deftypefun int gsl_linalg_hermtd_unpack (const gsl_matrix_complex * @var{A}, const gsl_vector_complex * @var{tau}, gsl_matrix_complex * @var{Q}, gsl_vector * @var{diag}, gsl_vector * @var{subdiag})
634 This function unpacks the encoded tridiagonal decomposition (@var{A},
635 @var{tau}) obtained from @code{gsl_linalg_hermtd_decomp} into the
636 unitary matrix @var{U}, the real vector of diagonal elements @var{diag} and
637 the real vector of subdiagonal elements @var{subdiag}. 
638 @end deftypefun
639
640 @deftypefun int gsl_linalg_hermtd_unpack_T (const gsl_matrix_complex * @var{A}, gsl_vector * @var{diag}, gsl_vector * @var{subdiag})
641 This function unpacks the diagonal and subdiagonal of the encoded
642 tridiagonal decomposition (@var{A}, @var{tau}) obtained from the
643 @code{gsl_linalg_hermtd_decomp} into the real vectors
644 @var{diag} and @var{subdiag}.
645 @end deftypefun
646
647 @node Hessenberg Decomposition of Real Matrices
648 @section Hessenberg Decomposition of Real Matrices
649 @cindex hessenberg decomposition
650
651 A general real matrix @math{A} can be decomposed by orthogonal
652 similarity transformations into the form
653 @tex
654 \beforedisplay
655 $$
656 A = U H U^T
657 $$
658 \afterdisplay
659 @end tex
660 @ifinfo
661
662 @example
663 A = U H U^T
664 @end example
665
666 @end ifinfo
667 where @math{U} is orthogonal and @math{H} is an upper Hessenberg matrix,
668 meaning that it has zeros below the first subdiagonal. The
669 Hessenberg reduction is the first step in the Schur decomposition
670 for the nonsymmetric eigenvalue problem, but has applications in
671 other areas as well.
672
673 @deftypefun int gsl_linalg_hessenberg_decomp (gsl_matrix * @var{A}, gsl_vector * @var{tau})
674 This function computes the Hessenberg decomposition of the matrix
675 @var{A} by applying the similarity transformation @math{H = U^T A U}.
676 On output, @math{H} is stored in the upper portion of @var{A}. The
677 information required to construct the matrix @math{U} is stored in
678 the lower triangular portion of @var{A}. @math{U} is a product
679 of @math{N - 2} Householder matrices. The Householder vectors
680 are stored in the lower portion of @var{A} (below the subdiagonal)
681 and the Householder coefficients are stored in the vector @var{tau}.
682 @var{tau} must be of length @var{N}.
683 @end deftypefun
684
685 @deftypefun int gsl_linalg_hessenberg_unpack (gsl_matrix * @var{H}, gsl_vector * @var{tau}, gsl_matrix * @var{U})
686 This function constructs the orthogonal matrix @math{U} from the
687 information stored in the Hessenberg matrix @var{H} along with the
688 vector @var{tau}. @var{H} and @var{tau} are outputs from
689 @code{gsl_linalg_hessenberg_decomp}.
690 @end deftypefun
691
692 @deftypefun int gsl_linalg_hessenberg_unpack_accum (gsl_matrix * @var{H}, gsl_vector * @var{tau}, gsl_matrix * @var{V})
693 This function is similar to @code{gsl_linalg_hessenberg_unpack}, except
694 it accumulates the matrix @var{U} into @var{V}, so that @math{V' = VU}.
695 The matrix @var{V} must be initialized prior to calling this function.
696 Setting @var{V} to the identity matrix provides the same result as
697 @code{gsl_linalg_hessenberg_unpack}. If @var{H} is order @var{N}, then
698 @var{V} must have @var{N} columns but may have any number of rows.
699 @end deftypefun
700
701 @deftypefun int gsl_linalg_hessenberg_set_zero (gsl_matrix * @var{H})
702 This function sets the lower triangular portion of @var{H}, below
703 the subdiagonal, to zero. It is useful for clearing out the
704 Householder vectors after calling @code{gsl_linalg_hessenberg_decomp}.
705 @end deftypefun
706
707 @node Hessenberg-Triangular Decomposition of Real Matrices
708 @section Hessenberg-Triangular Decomposition of Real Matrices
709 @cindex hessenberg triangular decomposition
710
711 A general real matrix pair (@math{A}, @math{B}) can be decomposed by
712 orthogonal similarity transformations into the form
713 @tex
714 \beforedisplay
715 $$
716 A = U H V^T
717 $$
718 $$
719 B = U R V^T
720 $$
721 \afterdisplay
722 @end tex
723 @ifinfo
724
725 @example
726 A = U H V^T
727 B = U R V^T
728 @end example
729
730 @end ifinfo
731 where @math{U} and @math{V} are orthogonal, @math{H} is an upper
732 Hessenberg matrix, and @math{R} is upper triangular. The
733 Hessenberg-Triangular reduction is the first step in the generalized
734 Schur decomposition for the generalized eigenvalue problem.
735
736 @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})
737 This function computes the Hessenberg-Triangular decomposition of the
738 matrix pair (@var{A}, @var{B}). On output, @math{H} is stored in @var{A},
739 and @math{R} is stored in @var{B}. If @var{U} and @var{V} are provided
740 (they may be null), the similarity transformations are stored in them.
741 Additional workspace of length @math{N} is needed in @var{work}.
742 @end deftypefun
743
744 @node Bidiagonalization
745 @section Bidiagonalization 
746 @cindex bidiagonalization of real matrices
747
748 A general matrix @math{A} can be factorized by similarity
749 transformations into the form,
750 @tex
751 \beforedisplay
752 $$
753 A = U B V^T
754 $$
755 \afterdisplay
756 @end tex
757 @ifinfo
758
759 @example
760 A = U B V^T
761 @end example
762
763 @end ifinfo
764 @noindent
765 where @math{U} and @math{V} are orthogonal matrices and @math{B} is a
766 @math{N}-by-@math{N} bidiagonal matrix with non-zero entries only on the
767 diagonal and superdiagonal.  The size of @var{U} is @math{M}-by-@math{N}
768 and the size of @var{V} is @math{N}-by-@math{N}.
769
770 @deftypefun int gsl_linalg_bidiag_decomp (gsl_matrix * @var{A}, gsl_vector * @var{tau_U}, gsl_vector * @var{tau_V})
771 This function factorizes the @math{M}-by-@math{N} matrix @var{A} into
772 bidiagonal form @math{U B V^T}.  The diagonal and superdiagonal of the
773 matrix @math{B} are stored in the diagonal and superdiagonal of @var{A}.
774 The orthogonal matrices @math{U} and @var{V} are stored as compressed
775 Householder vectors in the remaining elements of @var{A}.  The
776 Householder coefficients are stored in the vectors @var{tau_U} and
777 @var{tau_V}.  The length of @var{tau_U} must equal the number of
778 elements in the diagonal of @var{A} and the length of @var{tau_V} should
779 be one element shorter.
780 @end deftypefun
781
782 @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})
783 This function unpacks the bidiagonal decomposition of @var{A} given by
784 @code{gsl_linalg_bidiag_decomp}, (@var{A}, @var{tau_U}, @var{tau_V})
785 into the separate orthogonal matrices @var{U}, @var{V} and the diagonal
786 vector @var{diag} and superdiagonal @var{superdiag}.  Note that @var{U}
787 is stored as a compact @math{M}-by-@math{N} orthogonal matrix satisfying
788 @math{U^T U = I} for efficiency.
789 @end deftypefun
790
791 @deftypefun int gsl_linalg_bidiag_unpack2 (gsl_matrix * @var{A}, gsl_vector * @var{tau_U}, gsl_vector * @var{tau_V}, gsl_matrix * @var{V})
792 This function unpacks the bidiagonal decomposition of @var{A} given by
793 @code{gsl_linalg_bidiag_decomp}, (@var{A}, @var{tau_U}, @var{tau_V})
794 into the separate orthogonal matrices @var{U}, @var{V} and the diagonal
795 vector @var{diag} and superdiagonal @var{superdiag}.  The matrix @var{U}
796 is stored in-place in @var{A}.
797 @end deftypefun
798
799 @deftypefun int gsl_linalg_bidiag_unpack_B (const gsl_matrix * @var{A}, gsl_vector * @var{diag}, gsl_vector * @var{superdiag})
800 This function unpacks the diagonal and superdiagonal of the bidiagonal
801 decomposition of @var{A} given by @code{gsl_linalg_bidiag_decomp}, into
802 the diagonal vector @var{diag} and superdiagonal vector @var{superdiag}.
803 @end deftypefun
804
805 @node Householder Transformations
806 @section Householder Transformations
807 @cindex Householder matrix
808 @cindex Householder transformation
809 @cindex transformation, Householder
810
811 A Householder transformation is a rank-1 modification of the identity
812 matrix which can be used to zero out selected elements of a vector.  A
813 Householder matrix @math{P} takes the form,
814 @tex
815 \beforedisplay
816 $$
817 P = I - \tau v v^T
818 $$
819 \afterdisplay
820 @end tex
821 @ifinfo
822
823 @example
824 P = I - \tau v v^T
825 @end example
826
827 @end ifinfo
828 @noindent
829 where @math{v} is a vector (called the @dfn{Householder vector}) and
830 @math{\tau = 2/(v^T v)}.  The functions described in this section use the
831 rank-1 structure of the Householder matrix to create and apply
832 Householder transformations efficiently.
833
834 @deftypefun double gsl_linalg_householder_transform (gsl_vector * @var{v})
835 @deftypefunx gsl_complex gsl_linalg_complex_householder_transform (gsl_vector_complex * @var{v})
836 This function prepares a Householder transformation @math{P = I - \tau v
837 v^T} which can be used to zero all the elements of the input vector except
838 the first.  On output the transformation is stored in the vector @var{v}
839 and the scalar @math{\tau} is returned.
840 @end deftypefun
841
842 @deftypefun int gsl_linalg_householder_hm (double tau, const gsl_vector * v, gsl_matrix * A)
843 @deftypefunx int gsl_linalg_complex_householder_hm (gsl_complex tau, const gsl_vector_complex * v, gsl_matrix_complex * A)
844 This function applies the Householder matrix @math{P} defined by the
845 scalar @var{tau} and the vector @var{v} to the left-hand side of the
846 matrix @var{A}. On output the result @math{P A} is stored in @var{A}.
847 @end deftypefun
848
849 @deftypefun int gsl_linalg_householder_mh (double tau, const gsl_vector * v, gsl_matrix * A)
850 @deftypefunx int gsl_linalg_complex_householder_mh (gsl_complex tau, const gsl_vector_complex * v, gsl_matrix_complex * A)
851 This function applies the Householder matrix @math{P} defined by the
852 scalar @var{tau} and the vector @var{v} to the right-hand side of the
853 matrix @var{A}. On output the result @math{A P} is stored in @var{A}.
854 @end deftypefun
855
856 @deftypefun int gsl_linalg_householder_hv (double tau, const gsl_vector * v, gsl_vector * w)
857 @deftypefunx int gsl_linalg_complex_householder_hv (gsl_complex tau, const gsl_vector_complex * v, gsl_vector_complex * w)
858 This function applies the Householder transformation @math{P} defined by
859 the scalar @var{tau} and the vector @var{v} to the vector @var{w}.  On
860 output the result @math{P w} is stored in @var{w}.
861 @end deftypefun
862
863 @comment @deftypefun int gsl_linalg_householder_hm1 (double tau, gsl_matrix * A)
864 @comment This function applies the Householder transform, defined by the scalar
865 @comment @var{tau} and the vector @var{v}, to a matrix being build up from the
866 @comment identity matrix, using the first column of @var{A} as a householder vector.
867 @comment @end deftypefun
868
869 @node Householder solver for linear systems
870 @section Householder solver for linear systems
871 @cindex solution of linear system by Householder transformations
872 @cindex Householder linear solver
873
874 @deftypefun int gsl_linalg_HH_solve (gsl_matrix * @var{A}, const gsl_vector * @var{b}, gsl_vector * @var{x})
875 This function solves the system @math{A x = b} directly using
876 Householder transformations. On output the solution is stored in @var{x}
877 and @var{b} is not modified. The matrix @var{A} is destroyed by the
878 Householder transformations.
879 @end deftypefun
880
881 @deftypefun int gsl_linalg_HH_svx (gsl_matrix * @var{A}, gsl_vector * @var{x})
882 This function solves the system @math{A x = b} in-place using
883 Householder transformations.  On input @var{x} should contain the
884 right-hand side @math{b}, which is replaced by the solution on output.  The
885 matrix @var{A} is destroyed by the Householder transformations.
886 @end deftypefun
887
888 @node Tridiagonal Systems
889 @section Tridiagonal Systems
890 @cindex tridiagonal systems
891
892 The functions described in this section efficiently solve symmetric,
893 non-symmetric and cyclic tridiagonal systems with minimal storage.
894 Note that the current implementations of these functions use a variant
895 of Cholesky decomposition, so the tridiagonal matrix must be positive
896 definite.  For non-positive definite matrices, the functions return
897 the error code @code{GSL_ESING}.
898
899 @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})
900 This function solves the general @math{N}-by-@math{N} system @math{A x =
901 b} where @var{A} is tridiagonal (@c{$N\geq 2$}
902 @math{N >= 2}). The super-diagonal and
903 sub-diagonal vectors @var{e} and @var{f} must be one element shorter
904 than the diagonal vector @var{diag}.  The form of @var{A} for the 4-by-4
905 case is shown below,
906 @tex
907 \beforedisplay
908 $$
909 A = \pmatrix{d_0&e_0&  0& 0\cr
910              f_0&d_1&e_1& 0\cr
911              0  &f_1&d_2&e_2\cr 
912              0  &0  &f_2&d_3\cr}
913 $$
914 \afterdisplay
915 @end tex
916 @ifinfo
917
918 @example
919 A = ( d_0 e_0  0   0  )
920     ( f_0 d_1 e_1  0  )
921     (  0  f_1 d_2 e_2 )
922     (  0   0  f_2 d_3 )
923 @end example
924 @end ifinfo
925 @noindent
926 @end deftypefun
927
928 @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})
929 This function solves the general @math{N}-by-@math{N} system @math{A x =
930 b} where @var{A} is symmetric tridiagonal (@c{$N\geq 2$}
931 @math{N >= 2}).  The off-diagonal vector
932 @var{e} must be one element shorter than the diagonal vector @var{diag}.
933 The form of @var{A} for the 4-by-4 case is shown below,
934 @tex
935 \beforedisplay
936 $$
937 A = \pmatrix{d_0&e_0&  0& 0\cr
938              e_0&d_1&e_1& 0\cr
939              0  &e_1&d_2&e_2\cr 
940              0  &0  &e_2&d_3\cr}
941 $$
942 \afterdisplay
943 @end tex
944 @ifinfo
945
946 @example
947 A = ( d_0 e_0  0   0  )
948     ( e_0 d_1 e_1  0  )
949     (  0  e_1 d_2 e_2 )
950     (  0   0  e_2 d_3 )
951 @end example
952 @end ifinfo
953 @end deftypefun
954
955 @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})
956 This function solves the general @math{N}-by-@math{N} system @math{A x =
957 b} where @var{A} is cyclic tridiagonal (@c{$N\geq 3$}
958 @math{N >= 3}).  The cyclic super-diagonal and
959 sub-diagonal vectors @var{e} and @var{f} must have the same number of
960 elements as the diagonal vector @var{diag}.  The form of @var{A} for the
961 4-by-4 case is shown below,
962 @tex
963 \beforedisplay
964 $$
965 A = \pmatrix{d_0&e_0& 0 &f_3\cr
966              f_0&d_1&e_1& 0 \cr
967               0 &f_1&d_2&e_2\cr 
968              e_3& 0 &f_2&d_3\cr}
969 $$
970 \afterdisplay
971 @end tex
972 @ifinfo
973
974 @example
975 A = ( d_0 e_0  0  f_3 )
976     ( f_0 d_1 e_1  0  )
977     (  0  f_1 d_2 e_2 )
978     ( e_3  0  f_2 d_3 )
979 @end example
980 @end ifinfo
981 @end deftypefun
982
983
984 @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})
985 This function solves the general @math{N}-by-@math{N} system @math{A x =
986 b} where @var{A} is symmetric cyclic tridiagonal (@c{$N\geq 3$}
987 @math{N >= 3}).  The cyclic
988 off-diagonal vector @var{e} must have the same number of elements as the
989 diagonal vector @var{diag}.  The form of @var{A} for the 4-by-4 case is
990 shown below,
991 @tex
992 \beforedisplay
993 $$
994 A = \pmatrix{d_0&e_0& 0 &e_3\cr
995              e_0&d_1&e_1& 0 \cr
996               0 &e_1&d_2&e_2\cr 
997              e_3& 0 &e_2&d_3\cr}
998 $$
999 \afterdisplay
1000 @end tex
1001 @ifinfo
1002
1003 @example
1004 A = ( d_0 e_0  0  e_3 )
1005     ( e_0 d_1 e_1  0  )
1006     (  0  e_1 d_2 e_2 )
1007     ( e_3  0  e_2 d_3 )
1008 @end example
1009 @end ifinfo
1010 @end deftypefun
1011
1012 @node Balancing
1013 @section Balancing
1014 @cindex balancing matrices
1015
1016 The process of balancing a matrix applies similarity transformations
1017 to make the rows and columns have comparable norms. This is
1018 useful, for example, to reduce roundoff errors in the solution
1019 of eigenvalue problems. Balancing a matrix @math{A} consists
1020 of replacing @math{A} with a similar matrix
1021 @tex
1022 \beforedisplay
1023 $$
1024 A' = D^{-1} A D
1025 $$
1026 \afterdisplay
1027 @end tex
1028 @ifinfo
1029
1030 @example
1031 A' = D^(-1) A D
1032 @end example
1033
1034 @end ifinfo
1035 where @math{D} is a diagonal matrix whose entries are powers
1036 of the floating point radix.
1037
1038 @deftypefun int gsl_linalg_balance_matrix (gsl_matrix * @var{A}, gsl_vector * @var{D})
1039 This function replaces the matrix @var{A} with its balanced counterpart
1040 and stores the diagonal elements of the similarity transformation
1041 into the vector @var{D}.
1042 @end deftypefun
1043
1044 @node Linear Algebra Examples
1045 @section Examples
1046
1047 The following program solves the linear system @math{A x = b}. The
1048 system to be solved is,
1049 @tex
1050 \beforedisplay
1051 $$
1052 \left(
1053 \matrix{0.18& 0.60& 0.57& 0.96\cr
1054 0.41& 0.24& 0.99& 0.58\cr
1055 0.14& 0.30& 0.97& 0.66\cr
1056 0.51& 0.13& 0.19& 0.85}
1057 \right)
1058 \left(
1059 \matrix{x_0\cr
1060 x_1\cr
1061 x_2\cr
1062 x_3}
1063 \right)
1064 =
1065 \left(
1066 \matrix{1.0\cr
1067 2.0\cr
1068 3.0\cr
1069 4.0}
1070 \right)
1071 $$
1072 \afterdisplay
1073 @end tex
1074 @ifinfo
1075
1076 @example
1077 [ 0.18 0.60 0.57 0.96 ] [x0]   [1.0]
1078 [ 0.41 0.24 0.99 0.58 ] [x1] = [2.0]
1079 [ 0.14 0.30 0.97 0.66 ] [x2]   [3.0]
1080 [ 0.51 0.13 0.19 0.85 ] [x3]   [4.0]
1081 @end example
1082
1083 @end ifinfo
1084 @noindent
1085 and the solution is found using LU decomposition of the matrix @math{A}.
1086
1087 @example
1088 @verbatiminclude examples/linalglu.c
1089 @end example
1090
1091 @noindent
1092 Here is the output from the program,
1093
1094 @example
1095 @verbatiminclude examples/linalglu.out
1096 @end example
1097
1098 @noindent
1099 This can be verified by multiplying the solution @math{x} by the
1100 original matrix @math{A} using @sc{gnu octave},
1101
1102 @example
1103 octave> A = [ 0.18, 0.60, 0.57, 0.96;
1104               0.41, 0.24, 0.99, 0.58; 
1105               0.14, 0.30, 0.97, 0.66; 
1106               0.51, 0.13, 0.19, 0.85 ];
1107
1108 octave> x = [ -4.05205; -12.6056; 1.66091; 8.69377];
1109
1110 octave> A * x
1111 ans =
1112   1.0000
1113   2.0000
1114   3.0000
1115   4.0000
1116 @end example
1117
1118 @noindent
1119 This reproduces the original right-hand side vector, @math{b}, in
1120 accordance with the equation @math{A x = b}.
1121
1122 @node Linear Algebra References and Further Reading
1123 @section References and Further Reading
1124
1125 Further information on the algorithms described in this section can be
1126 found in the following book,
1127
1128 @itemize @asis
1129 @item
1130 G. H. Golub, C. F. Van Loan, @cite{Matrix Computations} (3rd Ed, 1996),
1131 Johns Hopkins University Press, ISBN 0-8018-5414-8.
1132 @end itemize
1133
1134 @noindent
1135 The @sc{lapack} library is described in the following manual,
1136
1137 @itemize @asis
1138 @item
1139 @cite{LAPACK Users' Guide} (Third Edition, 1999), Published by SIAM,
1140 ISBN 0-89871-447-8.
1141
1142 @uref{http://www.netlib.org/lapack} 
1143 @end itemize
1144
1145 @noindent
1146 The @sc{lapack} source code can be found at the website above, along
1147 with an online copy of the users guide.
1148
1149 @noindent
1150 The Modified Golub-Reinsch algorithm is described in the following paper,
1151
1152 @itemize @asis
1153 @item
1154 T.F. Chan, ``An Improved Algorithm for Computing the Singular Value
1155 Decomposition'', @cite{ACM Transactions on Mathematical Software}, 8
1156 (1982), pp 72--83.
1157 @end itemize
1158
1159 @noindent
1160 The Jacobi algorithm for singular value decomposition is described in
1161 the following papers,
1162
1163 @itemize @asis
1164 @item
1165 J.C. Nash, ``A one-sided transformation method for the singular value
1166 decomposition and algebraic eigenproblem'', @cite{Computer Journal},
1167 Volume 18, Number 1 (1973), p 74--76
1168
1169 @item
1170 James Demmel, Kresimir Veselic, ``Jacobi's Method is more accurate than
1171 QR'', @cite{Lapack Working Note 15} (LAWN-15), October 1989. Available
1172 from netlib, @uref{http://www.netlib.org/lapack/} in the @code{lawns} or
1173 @code{lawnspdf} directories.
1174 @end itemize
1175
1176
1177