gsl_linalg__qr.c 14 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601
  1. /* linalg/qr.c
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 Gerard Jungman, Brian Gough
  4. *
  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.
  9. *
  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.
  14. *
  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.
  18. */
  19. /* Author: G. Jungman */
  20. #include "gsl__config.h"
  21. #include <stdlib.h>
  22. #include <string.h>
  23. #include "gsl_math.h"
  24. #include "gsl_vector.h"
  25. #include "gsl_matrix.h"
  26. #include "gsl_blas.h"
  27. #include "gsl_linalg.h"
  28. #define REAL double
  29. #include "gsl_linalg__givens.c"
  30. #include "gsl_linalg__apply_givens.c"
  31. /* Factorise a general M x N matrix A into
  32. *
  33. * A = Q R
  34. *
  35. * where Q is orthogonal (M x M) and R is upper triangular (M x N).
  36. *
  37. * Q is stored as a packed set of Householder transformations in the
  38. * strict lower triangular part of the input matrix.
  39. *
  40. * R is stored in the diagonal and upper triangle of the input matrix.
  41. *
  42. * The full matrix for Q can be obtained as the product
  43. *
  44. * Q = Q_k .. Q_2 Q_1
  45. *
  46. * where k = MIN(M,N) and
  47. *
  48. * Q_i = (I - tau_i * v_i * v_i')
  49. *
  50. * and where v_i is a Householder vector
  51. *
  52. * v_i = [1, m(i+1,i), m(i+2,i), ... , m(M,i)]
  53. *
  54. * This storage scheme is the same as in LAPACK. */
  55. int
  56. gsl_linalg_QR_decomp (gsl_matrix * A, gsl_vector * tau)
  57. {
  58. const size_t M = A->size1;
  59. const size_t N = A->size2;
  60. if (tau->size != GSL_MIN (M, N))
  61. {
  62. GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
  63. }
  64. else
  65. {
  66. size_t i;
  67. for (i = 0; i < GSL_MIN (M, N); i++)
  68. {
  69. /* Compute the Householder transformation to reduce the j-th
  70. column of the matrix to a multiple of the j-th unit vector */
  71. gsl_vector_view c_full = gsl_matrix_column (A, i);
  72. gsl_vector_view c = gsl_vector_subvector (&(c_full.vector), i, M-i);
  73. double tau_i = gsl_linalg_householder_transform (&(c.vector));
  74. gsl_vector_set (tau, i, tau_i);
  75. /* Apply the transformation to the remaining columns and
  76. update the norms */
  77. if (i + 1 < N)
  78. {
  79. gsl_matrix_view m = gsl_matrix_submatrix (A, i, i + 1, M - i, N - (i + 1));
  80. gsl_linalg_householder_hm (tau_i, &(c.vector), &(m.matrix));
  81. }
  82. }
  83. return GSL_SUCCESS;
  84. }
  85. }
  86. /* Solves the system A x = b using the QR factorisation,
  87. * R x = Q^T b
  88. *
  89. * to obtain x. Based on SLATEC code.
  90. */
  91. int
  92. gsl_linalg_QR_solve (const gsl_matrix * QR, const gsl_vector * tau, const gsl_vector * b, gsl_vector * x)
  93. {
  94. if (QR->size1 != QR->size2)
  95. {
  96. GSL_ERROR ("QR matrix must be square", GSL_ENOTSQR);
  97. }
  98. else if (QR->size1 != b->size)
  99. {
  100. GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
  101. }
  102. else if (QR->size2 != x->size)
  103. {
  104. GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
  105. }
  106. else
  107. {
  108. /* Copy x <- b */
  109. gsl_vector_memcpy (x, b);
  110. /* Solve for x */
  111. gsl_linalg_QR_svx (QR, tau, x);
  112. return GSL_SUCCESS;
  113. }
  114. }
  115. /* Solves the system A x = b in place using the QR factorisation,
  116. * R x = Q^T b
  117. *
  118. * to obtain x. Based on SLATEC code.
  119. */
  120. int
  121. gsl_linalg_QR_svx (const gsl_matrix * QR, const gsl_vector * tau, gsl_vector * x)
  122. {
  123. if (QR->size1 != QR->size2)
  124. {
  125. GSL_ERROR ("QR matrix must be square", GSL_ENOTSQR);
  126. }
  127. else if (QR->size1 != x->size)
  128. {
  129. GSL_ERROR ("matrix size must match x/rhs size", GSL_EBADLEN);
  130. }
  131. else
  132. {
  133. /* compute rhs = Q^T b */
  134. gsl_linalg_QR_QTvec (QR, tau, x);
  135. /* Solve R x = rhs, storing x in-place */
  136. gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, QR, x);
  137. return GSL_SUCCESS;
  138. }
  139. }
  140. /* Find the least squares solution to the overdetermined system
  141. *
  142. * A x = b
  143. *
  144. * for M >= N using the QR factorization A = Q R.
  145. */
  146. int
  147. gsl_linalg_QR_lssolve (const gsl_matrix * QR, const gsl_vector * tau, const gsl_vector * b, gsl_vector * x, gsl_vector * residual)
  148. {
  149. const size_t M = QR->size1;
  150. const size_t N = QR->size2;
  151. if (M < N)
  152. {
  153. GSL_ERROR ("QR matrix must have M>=N", GSL_EBADLEN);
  154. }
  155. else if (M != b->size)
  156. {
  157. GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
  158. }
  159. else if (N != x->size)
  160. {
  161. GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
  162. }
  163. else if (M != residual->size)
  164. {
  165. GSL_ERROR ("matrix size must match residual size", GSL_EBADLEN);
  166. }
  167. else
  168. {
  169. gsl_matrix_const_view R = gsl_matrix_const_submatrix (QR, 0, 0, N, N);
  170. gsl_vector_view c = gsl_vector_subvector(residual, 0, N);
  171. gsl_vector_memcpy(residual, b);
  172. /* compute rhs = Q^T b */
  173. gsl_linalg_QR_QTvec (QR, tau, residual);
  174. /* Solve R x = rhs */
  175. gsl_vector_memcpy(x, &(c.vector));
  176. gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, &(R.matrix), x);
  177. /* Compute residual = b - A x = Q (Q^T b - R x) */
  178. gsl_vector_set_zero(&(c.vector));
  179. gsl_linalg_QR_Qvec(QR, tau, residual);
  180. return GSL_SUCCESS;
  181. }
  182. }
  183. int
  184. gsl_linalg_QR_Rsolve (const gsl_matrix * QR, const gsl_vector * b, gsl_vector * x)
  185. {
  186. if (QR->size1 != QR->size2)
  187. {
  188. GSL_ERROR ("QR matrix must be square", GSL_ENOTSQR);
  189. }
  190. else if (QR->size1 != b->size)
  191. {
  192. GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
  193. }
  194. else if (QR->size2 != x->size)
  195. {
  196. GSL_ERROR ("matrix size must match x size", GSL_EBADLEN);
  197. }
  198. else
  199. {
  200. /* Copy x <- b */
  201. gsl_vector_memcpy (x, b);
  202. /* Solve R x = b, storing x in-place */
  203. gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, QR, x);
  204. return GSL_SUCCESS;
  205. }
  206. }
  207. int
  208. gsl_linalg_QR_Rsvx (const gsl_matrix * QR, gsl_vector * x)
  209. {
  210. if (QR->size1 != QR->size2)
  211. {
  212. GSL_ERROR ("QR matrix must be square", GSL_ENOTSQR);
  213. }
  214. else if (QR->size1 != x->size)
  215. {
  216. GSL_ERROR ("matrix size must match rhs size", GSL_EBADLEN);
  217. }
  218. else
  219. {
  220. /* Solve R x = b, storing x in-place */
  221. gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, QR, x);
  222. return GSL_SUCCESS;
  223. }
  224. }
  225. int
  226. gsl_linalg_R_solve (const gsl_matrix * R, const gsl_vector * b, gsl_vector * x)
  227. {
  228. if (R->size1 != R->size2)
  229. {
  230. GSL_ERROR ("R matrix must be square", GSL_ENOTSQR);
  231. }
  232. else if (R->size1 != b->size)
  233. {
  234. GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
  235. }
  236. else if (R->size2 != x->size)
  237. {
  238. GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
  239. }
  240. else
  241. {
  242. /* Copy x <- b */
  243. gsl_vector_memcpy (x, b);
  244. /* Solve R x = b, storing x inplace in b */
  245. gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, R, x);
  246. return GSL_SUCCESS;
  247. }
  248. }
  249. int
  250. gsl_linalg_R_svx (const gsl_matrix * R, gsl_vector * x)
  251. {
  252. if (R->size1 != R->size2)
  253. {
  254. GSL_ERROR ("R matrix must be square", GSL_ENOTSQR);
  255. }
  256. else if (R->size2 != x->size)
  257. {
  258. GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
  259. }
  260. else
  261. {
  262. /* Solve R x = b, storing x inplace in b */
  263. gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, R, x);
  264. return GSL_SUCCESS;
  265. }
  266. }
  267. /* Form the product Q^T v from a QR factorized matrix
  268. */
  269. int
  270. gsl_linalg_QR_QTvec (const gsl_matrix * QR, const gsl_vector * tau, gsl_vector * v)
  271. {
  272. const size_t M = QR->size1;
  273. const size_t N = QR->size2;
  274. if (tau->size != GSL_MIN (M, N))
  275. {
  276. GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
  277. }
  278. else if (v->size != M)
  279. {
  280. GSL_ERROR ("vector size must be N", GSL_EBADLEN);
  281. }
  282. else
  283. {
  284. size_t i;
  285. /* compute Q^T v */
  286. for (i = 0; i < GSL_MIN (M, N); i++)
  287. {
  288. gsl_vector_const_view c = gsl_matrix_const_column (QR, i);
  289. gsl_vector_const_view h = gsl_vector_const_subvector (&(c.vector), i, M - i);
  290. gsl_vector_view w = gsl_vector_subvector (v, i, M - i);
  291. double ti = gsl_vector_get (tau, i);
  292. gsl_linalg_householder_hv (ti, &(h.vector), &(w.vector));
  293. }
  294. return GSL_SUCCESS;
  295. }
  296. }
  297. int
  298. gsl_linalg_QR_Qvec (const gsl_matrix * QR, const gsl_vector * tau, gsl_vector * v)
  299. {
  300. const size_t M = QR->size1;
  301. const size_t N = QR->size2;
  302. if (tau->size != GSL_MIN (M, N))
  303. {
  304. GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
  305. }
  306. else if (v->size != M)
  307. {
  308. GSL_ERROR ("vector size must be N", GSL_EBADLEN);
  309. }
  310. else
  311. {
  312. size_t i;
  313. /* compute Q^T v */
  314. for (i = GSL_MIN (M, N); i > 0 && i--;)
  315. {
  316. gsl_vector_const_view c = gsl_matrix_const_column (QR, i);
  317. gsl_vector_const_view h = gsl_vector_const_subvector (&(c.vector),
  318. i, M - i);
  319. gsl_vector_view w = gsl_vector_subvector (v, i, M - i);
  320. double ti = gsl_vector_get (tau, i);
  321. gsl_linalg_householder_hv (ti, &h.vector, &w.vector);
  322. }
  323. return GSL_SUCCESS;
  324. }
  325. }
  326. /* Form the product Q^T A from a QR factorized matrix */
  327. int
  328. gsl_linalg_QR_QTmat (const gsl_matrix * QR, const gsl_vector * tau, gsl_matrix * A)
  329. {
  330. const size_t M = QR->size1;
  331. const size_t N = QR->size2;
  332. if (tau->size != GSL_MIN (M, N))
  333. {
  334. GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
  335. }
  336. else if (A->size1 != M)
  337. {
  338. GSL_ERROR ("matrix must have M rows", GSL_EBADLEN);
  339. }
  340. else
  341. {
  342. size_t i;
  343. /* compute Q^T A */
  344. for (i = 0; i < GSL_MIN (M, N); i++)
  345. {
  346. gsl_vector_const_view c = gsl_matrix_const_column (QR, i);
  347. gsl_vector_const_view h = gsl_vector_const_subvector (&(c.vector), i, M - i);
  348. gsl_matrix_view m = gsl_matrix_submatrix(A, i, 0, M - i, A->size2);
  349. double ti = gsl_vector_get (tau, i);
  350. gsl_linalg_householder_hm (ti, &(h.vector), &(m.matrix));
  351. }
  352. return GSL_SUCCESS;
  353. }
  354. }
  355. /* Form the orthogonal matrix Q from the packed QR matrix */
  356. int
  357. gsl_linalg_QR_unpack (const gsl_matrix * QR, const gsl_vector * tau, gsl_matrix * Q, gsl_matrix * R)
  358. {
  359. const size_t M = QR->size1;
  360. const size_t N = QR->size2;
  361. if (Q->size1 != M || Q->size2 != M)
  362. {
  363. GSL_ERROR ("Q matrix must be M x M", GSL_ENOTSQR);
  364. }
  365. else if (R->size1 != M || R->size2 != N)
  366. {
  367. GSL_ERROR ("R matrix must be M x N", GSL_ENOTSQR);
  368. }
  369. else if (tau->size != GSL_MIN (M, N))
  370. {
  371. GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
  372. }
  373. else
  374. {
  375. size_t i, j;
  376. /* Initialize Q to the identity */
  377. gsl_matrix_set_identity (Q);
  378. for (i = GSL_MIN (M, N); i > 0 && i--;)
  379. {
  380. gsl_vector_const_view c = gsl_matrix_const_column (QR, i);
  381. gsl_vector_const_view h = gsl_vector_const_subvector (&c.vector,
  382. i, M - i);
  383. gsl_matrix_view m = gsl_matrix_submatrix (Q, i, i, M - i, M - i);
  384. double ti = gsl_vector_get (tau, i);
  385. gsl_linalg_householder_hm (ti, &h.vector, &m.matrix);
  386. }
  387. /* Form the right triangular matrix R from a packed QR matrix */
  388. for (i = 0; i < M; i++)
  389. {
  390. for (j = 0; j < i && j < N; j++)
  391. gsl_matrix_set (R, i, j, 0.0);
  392. for (j = i; j < N; j++)
  393. gsl_matrix_set (R, i, j, gsl_matrix_get (QR, i, j));
  394. }
  395. return GSL_SUCCESS;
  396. }
  397. }
  398. /* Update a QR factorisation for A= Q R , A' = A + u v^T,
  399. * Q' R' = QR + u v^T
  400. * = Q (R + Q^T u v^T)
  401. * = Q (R + w v^T)
  402. *
  403. * where w = Q^T u.
  404. *
  405. * Algorithm from Golub and Van Loan, "Matrix Computations", Section
  406. * 12.5 (Updating Matrix Factorizations, Rank-One Changes)
  407. */
  408. int
  409. gsl_linalg_QR_update (gsl_matrix * Q, gsl_matrix * R,
  410. gsl_vector * w, const gsl_vector * v)
  411. {
  412. const size_t M = R->size1;
  413. const size_t N = R->size2;
  414. if (Q->size1 != M || Q->size2 != M)
  415. {
  416. GSL_ERROR ("Q matrix must be M x M if R is M x N", GSL_ENOTSQR);
  417. }
  418. else if (w->size != M)
  419. {
  420. GSL_ERROR ("w must be length M if R is M x N", GSL_EBADLEN);
  421. }
  422. else if (v->size != N)
  423. {
  424. GSL_ERROR ("v must be length N if R is M x N", GSL_EBADLEN);
  425. }
  426. else
  427. {
  428. size_t j, k;
  429. double w0;
  430. /* Apply Given's rotations to reduce w to (|w|, 0, 0, ... , 0)
  431. J_1^T .... J_(n-1)^T w = +/- |w| e_1
  432. simultaneously applied to R, H = J_1^T ... J^T_(n-1) R
  433. so that H is upper Hessenberg. (12.5.2) */
  434. for (k = M - 1; k > 0; k--)
  435. {
  436. double c, s;
  437. double wk = gsl_vector_get (w, k);
  438. double wkm1 = gsl_vector_get (w, k - 1);
  439. create_givens (wkm1, wk, &c, &s);
  440. apply_givens_vec (w, k - 1, k, c, s);
  441. apply_givens_qr (M, N, Q, R, k - 1, k, c, s);
  442. }
  443. w0 = gsl_vector_get (w, 0);
  444. /* Add in w v^T (Equation 12.5.3) */
  445. for (j = 0; j < N; j++)
  446. {
  447. double r0j = gsl_matrix_get (R, 0, j);
  448. double vj = gsl_vector_get (v, j);
  449. gsl_matrix_set (R, 0, j, r0j + w0 * vj);
  450. }
  451. /* Apply Givens transformations R' = G_(n-1)^T ... G_1^T H
  452. Equation 12.5.4 */
  453. for (k = 1; k < GSL_MIN(M,N+1); k++)
  454. {
  455. double c, s;
  456. double diag = gsl_matrix_get (R, k - 1, k - 1);
  457. double offdiag = gsl_matrix_get (R, k, k - 1);
  458. create_givens (diag, offdiag, &c, &s);
  459. apply_givens_qr (M, N, Q, R, k - 1, k, c, s);
  460. gsl_matrix_set (R, k, k - 1, 0.0); /* exact zero of G^T */
  461. }
  462. return GSL_SUCCESS;
  463. }
  464. }
  465. int
  466. gsl_linalg_QR_QRsolve (gsl_matrix * Q, gsl_matrix * R, const gsl_vector * b, gsl_vector * x)
  467. {
  468. const size_t M = R->size1;
  469. const size_t N = R->size2;
  470. if (M != N)
  471. {
  472. return GSL_ENOTSQR;
  473. }
  474. else if (Q->size1 != M || b->size != M || x->size != M)
  475. {
  476. return GSL_EBADLEN;
  477. }
  478. else
  479. {
  480. /* compute sol = Q^T b */
  481. gsl_blas_dgemv (CblasTrans, 1.0, Q, b, 0.0, x);
  482. /* Solve R x = sol, storing x in-place */
  483. gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, R, x);
  484. return GSL_SUCCESS;
  485. }
  486. }