3 * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 Gerard Jungman, Brian Gough
5 * This program is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 3 of the License, or (at
8 * your option) any later version.
10 * This program is distributed in the hope that it will be useful, but
11 * WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13 * General Public License for more details.
15 * You should have received a copy of the GNU General Public License
16 * along with this program; if not, write to the Free Software
17 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
20 /* Author: G. Jungman */
25 #include <gsl/gsl_math.h>
26 #include <gsl/gsl_vector.h>
27 #include <gsl/gsl_matrix.h>
28 #include <gsl/gsl_blas.h>
30 #include <gsl/gsl_linalg.h>
35 #include "apply_givens.c"
37 /* Factorise a general M x N matrix A into
41 * where Q is orthogonal (M x M) and R is upper triangular (M x N).
43 * Q is stored as a packed set of Householder transformations in the
44 * strict lower triangular part of the input matrix.
46 * R is stored in the diagonal and upper triangle of the input matrix.
48 * The full matrix for Q can be obtained as the product
52 * where k = MIN(M,N) and
54 * Q_i = (I - tau_i * v_i * v_i')
56 * and where v_i is a Householder vector
58 * v_i = [1, m(i+1,i), m(i+2,i), ... , m(M,i)]
60 * This storage scheme is the same as in LAPACK. */
63 gsl_linalg_QR_decomp (gsl_matrix * A, gsl_vector * tau)
65 const size_t M = A->size1;
66 const size_t N = A->size2;
68 if (tau->size != GSL_MIN (M, N))
70 GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
76 for (i = 0; i < GSL_MIN (M, N); i++)
78 /* Compute the Householder transformation to reduce the j-th
79 column of the matrix to a multiple of the j-th unit vector */
81 gsl_vector_view c_full = gsl_matrix_column (A, i);
82 gsl_vector_view c = gsl_vector_subvector (&(c_full.vector), i, M-i);
84 double tau_i = gsl_linalg_householder_transform (&(c.vector));
86 gsl_vector_set (tau, i, tau_i);
88 /* Apply the transformation to the remaining columns and
93 gsl_matrix_view m = gsl_matrix_submatrix (A, i, i + 1, M - i, N - (i + 1));
94 gsl_linalg_householder_hm (tau_i, &(c.vector), &(m.matrix));
102 /* Solves the system A x = b using the QR factorisation,
106 * to obtain x. Based on SLATEC code.
110 gsl_linalg_QR_solve (const gsl_matrix * QR, const gsl_vector * tau, const gsl_vector * b, gsl_vector * x)
112 if (QR->size1 != QR->size2)
114 GSL_ERROR ("QR matrix must be square", GSL_ENOTSQR);
116 else if (QR->size1 != b->size)
118 GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
120 else if (QR->size2 != x->size)
122 GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
128 gsl_vector_memcpy (x, b);
132 gsl_linalg_QR_svx (QR, tau, x);
138 /* Solves the system A x = b in place using the QR factorisation,
142 * to obtain x. Based on SLATEC code.
146 gsl_linalg_QR_svx (const gsl_matrix * QR, const gsl_vector * tau, gsl_vector * x)
149 if (QR->size1 != QR->size2)
151 GSL_ERROR ("QR matrix must be square", GSL_ENOTSQR);
153 else if (QR->size1 != x->size)
155 GSL_ERROR ("matrix size must match x/rhs size", GSL_EBADLEN);
159 /* compute rhs = Q^T b */
161 gsl_linalg_QR_QTvec (QR, tau, x);
163 /* Solve R x = rhs, storing x in-place */
165 gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, QR, x);
172 /* Find the least squares solution to the overdetermined system
176 * for M >= N using the QR factorization A = Q R.
180 gsl_linalg_QR_lssolve (const gsl_matrix * QR, const gsl_vector * tau, const gsl_vector * b, gsl_vector * x, gsl_vector * residual)
182 const size_t M = QR->size1;
183 const size_t N = QR->size2;
187 GSL_ERROR ("QR matrix must have M>=N", GSL_EBADLEN);
189 else if (M != b->size)
191 GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
193 else if (N != x->size)
195 GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
197 else if (M != residual->size)
199 GSL_ERROR ("matrix size must match residual size", GSL_EBADLEN);
203 gsl_matrix_const_view R = gsl_matrix_const_submatrix (QR, 0, 0, N, N);
204 gsl_vector_view c = gsl_vector_subvector(residual, 0, N);
206 gsl_vector_memcpy(residual, b);
208 /* compute rhs = Q^T b */
210 gsl_linalg_QR_QTvec (QR, tau, residual);
212 /* Solve R x = rhs */
214 gsl_vector_memcpy(x, &(c.vector));
216 gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, &(R.matrix), x);
218 /* Compute residual = b - A x = Q (Q^T b - R x) */
220 gsl_vector_set_zero(&(c.vector));
222 gsl_linalg_QR_Qvec(QR, tau, residual);
230 gsl_linalg_QR_Rsolve (const gsl_matrix * QR, const gsl_vector * b, gsl_vector * x)
232 if (QR->size1 != QR->size2)
234 GSL_ERROR ("QR matrix must be square", GSL_ENOTSQR);
236 else if (QR->size1 != b->size)
238 GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
240 else if (QR->size2 != x->size)
242 GSL_ERROR ("matrix size must match x size", GSL_EBADLEN);
248 gsl_vector_memcpy (x, b);
250 /* Solve R x = b, storing x in-place */
252 gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, QR, x);
260 gsl_linalg_QR_Rsvx (const gsl_matrix * QR, gsl_vector * x)
262 if (QR->size1 != QR->size2)
264 GSL_ERROR ("QR matrix must be square", GSL_ENOTSQR);
266 else if (QR->size1 != x->size)
268 GSL_ERROR ("matrix size must match rhs size", GSL_EBADLEN);
272 /* Solve R x = b, storing x in-place */
274 gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, QR, x);
281 gsl_linalg_R_solve (const gsl_matrix * R, const gsl_vector * b, gsl_vector * x)
283 if (R->size1 != R->size2)
285 GSL_ERROR ("R matrix must be square", GSL_ENOTSQR);
287 else if (R->size1 != b->size)
289 GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
291 else if (R->size2 != x->size)
293 GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
299 gsl_vector_memcpy (x, b);
301 /* Solve R x = b, storing x inplace in b */
303 gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, R, x);
310 gsl_linalg_R_svx (const gsl_matrix * R, gsl_vector * x)
312 if (R->size1 != R->size2)
314 GSL_ERROR ("R matrix must be square", GSL_ENOTSQR);
316 else if (R->size2 != x->size)
318 GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
322 /* Solve R x = b, storing x inplace in b */
324 gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, R, x);
332 /* Form the product Q^T v from a QR factorized matrix
336 gsl_linalg_QR_QTvec (const gsl_matrix * QR, const gsl_vector * tau, gsl_vector * v)
338 const size_t M = QR->size1;
339 const size_t N = QR->size2;
341 if (tau->size != GSL_MIN (M, N))
343 GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
345 else if (v->size != M)
347 GSL_ERROR ("vector size must be N", GSL_EBADLEN);
355 for (i = 0; i < GSL_MIN (M, N); i++)
357 gsl_vector_const_view c = gsl_matrix_const_column (QR, i);
358 gsl_vector_const_view h = gsl_vector_const_subvector (&(c.vector), i, M - i);
359 gsl_vector_view w = gsl_vector_subvector (v, i, M - i);
360 double ti = gsl_vector_get (tau, i);
361 gsl_linalg_householder_hv (ti, &(h.vector), &(w.vector));
369 gsl_linalg_QR_Qvec (const gsl_matrix * QR, const gsl_vector * tau, gsl_vector * v)
371 const size_t M = QR->size1;
372 const size_t N = QR->size2;
374 if (tau->size != GSL_MIN (M, N))
376 GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
378 else if (v->size != M)
380 GSL_ERROR ("vector size must be N", GSL_EBADLEN);
388 for (i = GSL_MIN (M, N); i > 0 && i--;)
390 gsl_vector_const_view c = gsl_matrix_const_column (QR, i);
391 gsl_vector_const_view h = gsl_vector_const_subvector (&(c.vector),
393 gsl_vector_view w = gsl_vector_subvector (v, i, M - i);
394 double ti = gsl_vector_get (tau, i);
395 gsl_linalg_householder_hv (ti, &h.vector, &w.vector);
401 /* Form the product Q^T A from a QR factorized matrix */
404 gsl_linalg_QR_QTmat (const gsl_matrix * QR, const gsl_vector * tau, gsl_matrix * A)
406 const size_t M = QR->size1;
407 const size_t N = QR->size2;
409 if (tau->size != GSL_MIN (M, N))
411 GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
413 else if (A->size1 != M)
415 GSL_ERROR ("matrix must have M rows", GSL_EBADLEN);
423 for (i = 0; i < GSL_MIN (M, N); i++)
425 gsl_vector_const_view c = gsl_matrix_const_column (QR, i);
426 gsl_vector_const_view h = gsl_vector_const_subvector (&(c.vector), i, M - i);
427 gsl_matrix_view m = gsl_matrix_submatrix(A, i, 0, M - i, A->size2);
428 double ti = gsl_vector_get (tau, i);
429 gsl_linalg_householder_hm (ti, &(h.vector), &(m.matrix));
436 /* Form the orthogonal matrix Q from the packed QR matrix */
439 gsl_linalg_QR_unpack (const gsl_matrix * QR, const gsl_vector * tau, gsl_matrix * Q, gsl_matrix * R)
441 const size_t M = QR->size1;
442 const size_t N = QR->size2;
444 if (Q->size1 != M || Q->size2 != M)
446 GSL_ERROR ("Q matrix must be M x M", GSL_ENOTSQR);
448 else if (R->size1 != M || R->size2 != N)
450 GSL_ERROR ("R matrix must be M x N", GSL_ENOTSQR);
452 else if (tau->size != GSL_MIN (M, N))
454 GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
460 /* Initialize Q to the identity */
462 gsl_matrix_set_identity (Q);
464 for (i = GSL_MIN (M, N); i > 0 && i--;)
466 gsl_vector_const_view c = gsl_matrix_const_column (QR, i);
467 gsl_vector_const_view h = gsl_vector_const_subvector (&c.vector,
469 gsl_matrix_view m = gsl_matrix_submatrix (Q, i, i, M - i, M - i);
470 double ti = gsl_vector_get (tau, i);
471 gsl_linalg_householder_hm (ti, &h.vector, &m.matrix);
474 /* Form the right triangular matrix R from a packed QR matrix */
476 for (i = 0; i < M; i++)
478 for (j = 0; j < i && j < N; j++)
479 gsl_matrix_set (R, i, j, 0.0);
481 for (j = i; j < N; j++)
482 gsl_matrix_set (R, i, j, gsl_matrix_get (QR, i, j));
490 /* Update a QR factorisation for A= Q R , A' = A + u v^T,
493 * = Q (R + Q^T u v^T)
498 * Algorithm from Golub and Van Loan, "Matrix Computations", Section
499 * 12.5 (Updating Matrix Factorizations, Rank-One Changes)
503 gsl_linalg_QR_update (gsl_matrix * Q, gsl_matrix * R,
504 gsl_vector * w, const gsl_vector * v)
506 const size_t M = R->size1;
507 const size_t N = R->size2;
509 if (Q->size1 != M || Q->size2 != M)
511 GSL_ERROR ("Q matrix must be M x M if R is M x N", GSL_ENOTSQR);
513 else if (w->size != M)
515 GSL_ERROR ("w must be length M if R is M x N", GSL_EBADLEN);
517 else if (v->size != N)
519 GSL_ERROR ("v must be length N if R is M x N", GSL_EBADLEN);
526 /* Apply Given's rotations to reduce w to (|w|, 0, 0, ... , 0)
528 J_1^T .... J_(n-1)^T w = +/- |w| e_1
530 simultaneously applied to R, H = J_1^T ... J^T_(n-1) R
531 so that H is upper Hessenberg. (12.5.2) */
533 for (k = M - 1; k > 0; k--)
536 double wk = gsl_vector_get (w, k);
537 double wkm1 = gsl_vector_get (w, k - 1);
539 create_givens (wkm1, wk, &c, &s);
540 apply_givens_vec (w, k - 1, k, c, s);
541 apply_givens_qr (M, N, Q, R, k - 1, k, c, s);
544 w0 = gsl_vector_get (w, 0);
546 /* Add in w v^T (Equation 12.5.3) */
548 for (j = 0; j < N; j++)
550 double r0j = gsl_matrix_get (R, 0, j);
551 double vj = gsl_vector_get (v, j);
552 gsl_matrix_set (R, 0, j, r0j + w0 * vj);
555 /* Apply Givens transformations R' = G_(n-1)^T ... G_1^T H
558 for (k = 1; k < GSL_MIN(M,N+1); k++)
561 double diag = gsl_matrix_get (R, k - 1, k - 1);
562 double offdiag = gsl_matrix_get (R, k, k - 1);
564 create_givens (diag, offdiag, &c, &s);
565 apply_givens_qr (M, N, Q, R, k - 1, k, c, s);
567 gsl_matrix_set (R, k, k - 1, 0.0); /* exact zero of G^T */
575 gsl_linalg_QR_QRsolve (gsl_matrix * Q, gsl_matrix * R, const gsl_vector * b, gsl_vector * x)
577 const size_t M = R->size1;
578 const size_t N = R->size2;
584 else if (Q->size1 != M || b->size != M || x->size != M)
590 /* compute sol = Q^T b */
592 gsl_blas_dgemv (CblasTrans, 1.0, Q, b, 0.0, x);
594 /* Solve R x = sol, storing x in-place */
596 gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, R, x);