gsl_linalg__qrpt.c 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487
  1. /* linalg/qrpt.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. #include "gsl__config.h"
  20. #include <stdlib.h>
  21. #include <string.h>
  22. #include "gsl_math.h"
  23. #include "gsl_vector.h"
  24. #include "gsl_matrix.h"
  25. #include "gsl_permute_vector.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 P = Q R
  34. *
  35. * where Q is orthogonal (M x M) and R is upper triangular (M x N).
  36. * When A is rank deficient, r = rank(A) < n, then the permutation is
  37. * used to ensure that the lower n - r rows of R are zero and the first
  38. * r columns of Q form an orthonormal basis for A.
  39. *
  40. * Q is stored as a packed set of Householder transformations in the
  41. * strict lower triangular part of the input matrix.
  42. *
  43. * R is stored in the diagonal and upper triangle of the input matrix.
  44. *
  45. * P: column j of P is column k of the identity matrix, where k =
  46. * permutation->data[j]
  47. *
  48. * The full matrix for Q can be obtained as the product
  49. *
  50. * Q = Q_k .. Q_2 Q_1
  51. *
  52. * where k = MIN(M,N) and
  53. *
  54. * Q_i = (I - tau_i * v_i * v_i')
  55. *
  56. * and where v_i is a Householder vector
  57. *
  58. * v_i = [1, m(i+1,i), m(i+2,i), ... , m(M,i)]
  59. *
  60. * This storage scheme is the same as in LAPACK. See LAPACK's
  61. * dgeqpf.f for details.
  62. *
  63. */
  64. int
  65. gsl_linalg_QRPT_decomp (gsl_matrix * A, gsl_vector * tau, gsl_permutation * p, int *signum, gsl_vector * norm)
  66. {
  67. const size_t M = A->size1;
  68. const size_t N = A->size2;
  69. if (tau->size != GSL_MIN (M, N))
  70. {
  71. GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
  72. }
  73. else if (p->size != N)
  74. {
  75. GSL_ERROR ("permutation size must be N", GSL_EBADLEN);
  76. }
  77. else if (norm->size != N)
  78. {
  79. GSL_ERROR ("norm size must be N", GSL_EBADLEN);
  80. }
  81. else
  82. {
  83. size_t i;
  84. *signum = 1;
  85. gsl_permutation_init (p); /* set to identity */
  86. /* Compute column norms and store in workspace */
  87. for (i = 0; i < N; i++)
  88. {
  89. gsl_vector_view c = gsl_matrix_column (A, i);
  90. double x = gsl_blas_dnrm2 (&c.vector);
  91. gsl_vector_set (norm, i, x);
  92. }
  93. for (i = 0; i < GSL_MIN (M, N); i++)
  94. {
  95. /* Bring the column of largest norm into the pivot position */
  96. double max_norm = gsl_vector_get(norm, i);
  97. size_t j, kmax = i;
  98. for (j = i + 1; j < N; j++)
  99. {
  100. double x = gsl_vector_get (norm, j);
  101. if (x > max_norm)
  102. {
  103. max_norm = x;
  104. kmax = j;
  105. }
  106. }
  107. if (kmax != i)
  108. {
  109. gsl_matrix_swap_columns (A, i, kmax);
  110. gsl_permutation_swap (p, i, kmax);
  111. gsl_vector_swap_elements(norm,i,kmax);
  112. (*signum) = -(*signum);
  113. }
  114. /* Compute the Householder transformation to reduce the j-th
  115. column of the matrix to a multiple of the j-th unit vector */
  116. {
  117. gsl_vector_view c_full = gsl_matrix_column (A, i);
  118. gsl_vector_view c = gsl_vector_subvector (&c_full.vector,
  119. i, M - i);
  120. double tau_i = gsl_linalg_householder_transform (&c.vector);
  121. gsl_vector_set (tau, i, tau_i);
  122. /* Apply the transformation to the remaining columns */
  123. if (i + 1 < N)
  124. {
  125. gsl_matrix_view m = gsl_matrix_submatrix (A, i, i + 1, M - i, N - (i+1));
  126. gsl_linalg_householder_hm (tau_i, &c.vector, &m.matrix);
  127. }
  128. }
  129. /* Update the norms of the remaining columns too */
  130. if (i + 1 < M)
  131. {
  132. for (j = i + 1; j < N; j++)
  133. {
  134. double x = gsl_vector_get (norm, j);
  135. if (x > 0.0)
  136. {
  137. double y = 0;
  138. double temp= gsl_matrix_get (A, i, j) / x;
  139. if (fabs (temp) >= 1)
  140. y = 0.0;
  141. else
  142. y = x * sqrt (1 - temp * temp);
  143. /* recompute norm to prevent loss of accuracy */
  144. if (fabs (y / x) < sqrt (20.0) * GSL_SQRT_DBL_EPSILON)
  145. {
  146. gsl_vector_view c_full = gsl_matrix_column (A, j);
  147. gsl_vector_view c =
  148. gsl_vector_subvector(&c_full.vector,
  149. i+1, M - (i+1));
  150. y = gsl_blas_dnrm2 (&c.vector);
  151. }
  152. gsl_vector_set (norm, j, y);
  153. }
  154. }
  155. }
  156. }
  157. return GSL_SUCCESS;
  158. }
  159. }
  160. int
  161. gsl_linalg_QRPT_decomp2 (const gsl_matrix * A, gsl_matrix * q, gsl_matrix * r, gsl_vector * tau, gsl_permutation * p, int *signum, gsl_vector * norm)
  162. {
  163. const size_t M = A->size1;
  164. const size_t N = A->size2;
  165. if (q->size1 != M || q->size2 !=M)
  166. {
  167. GSL_ERROR ("q must be M x M", GSL_EBADLEN);
  168. }
  169. else if (r->size1 != M || r->size2 !=N)
  170. {
  171. GSL_ERROR ("r must be M x N", GSL_EBADLEN);
  172. }
  173. else if (tau->size != GSL_MIN (M, N))
  174. {
  175. GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
  176. }
  177. else if (p->size != N)
  178. {
  179. GSL_ERROR ("permutation size must be N", GSL_EBADLEN);
  180. }
  181. else if (norm->size != N)
  182. {
  183. GSL_ERROR ("norm size must be N", GSL_EBADLEN);
  184. }
  185. gsl_matrix_memcpy (r, A);
  186. gsl_linalg_QRPT_decomp (r, tau, p, signum, norm);
  187. /* FIXME: aliased arguments depends on behavior of unpack routine! */
  188. gsl_linalg_QR_unpack (r, tau, q, r);
  189. return GSL_SUCCESS;
  190. }
  191. /* Solves the system A x = b using the Q R P^T factorisation,
  192. R z = Q^T b
  193. x = P z;
  194. to obtain x. Based on SLATEC code. */
  195. int
  196. gsl_linalg_QRPT_solve (const gsl_matrix * QR,
  197. const gsl_vector * tau,
  198. const gsl_permutation * p,
  199. const gsl_vector * b,
  200. gsl_vector * x)
  201. {
  202. if (QR->size1 != QR->size2)
  203. {
  204. GSL_ERROR ("QR matrix must be square", GSL_ENOTSQR);
  205. }
  206. else if (QR->size1 != p->size)
  207. {
  208. GSL_ERROR ("matrix size must match permutation size", GSL_EBADLEN);
  209. }
  210. else if (QR->size1 != b->size)
  211. {
  212. GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
  213. }
  214. else if (QR->size2 != x->size)
  215. {
  216. GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
  217. }
  218. else
  219. {
  220. gsl_vector_memcpy (x, b);
  221. gsl_linalg_QRPT_svx (QR, tau, p, x);
  222. return GSL_SUCCESS;
  223. }
  224. }
  225. int
  226. gsl_linalg_QRPT_svx (const gsl_matrix * QR,
  227. const gsl_vector * tau,
  228. const gsl_permutation * p,
  229. gsl_vector * x)
  230. {
  231. if (QR->size1 != QR->size2)
  232. {
  233. GSL_ERROR ("QR matrix must be square", GSL_ENOTSQR);
  234. }
  235. else if (QR->size1 != p->size)
  236. {
  237. GSL_ERROR ("matrix size must match permutation size", GSL_EBADLEN);
  238. }
  239. else if (QR->size2 != x->size)
  240. {
  241. GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
  242. }
  243. else
  244. {
  245. /* compute sol = Q^T b */
  246. gsl_linalg_QR_QTvec (QR, tau, x);
  247. /* Solve R x = sol, storing x inplace in sol */
  248. gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, QR, x);
  249. gsl_permute_vector_inverse (p, x);
  250. return GSL_SUCCESS;
  251. }
  252. }
  253. int
  254. gsl_linalg_QRPT_QRsolve (const gsl_matrix * Q, const gsl_matrix * R,
  255. const gsl_permutation * p,
  256. const gsl_vector * b,
  257. gsl_vector * x)
  258. {
  259. if (Q->size1 != Q->size2 || R->size1 != R->size2)
  260. {
  261. return GSL_ENOTSQR;
  262. }
  263. else if (Q->size1 != p->size || Q->size1 != R->size1
  264. || Q->size1 != b->size)
  265. {
  266. return GSL_EBADLEN;
  267. }
  268. else
  269. {
  270. /* compute b' = Q^T b */
  271. gsl_blas_dgemv (CblasTrans, 1.0, Q, b, 0.0, x);
  272. /* Solve R x = b', storing x inplace */
  273. gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, R, x);
  274. /* Apply permutation to solution in place */
  275. gsl_permute_vector_inverse (p, x);
  276. return GSL_SUCCESS;
  277. }
  278. }
  279. int
  280. gsl_linalg_QRPT_Rsolve (const gsl_matrix * QR,
  281. const gsl_permutation * p,
  282. const gsl_vector * b,
  283. gsl_vector * x)
  284. {
  285. if (QR->size1 != QR->size2)
  286. {
  287. GSL_ERROR ("QR matrix must be square", GSL_ENOTSQR);
  288. }
  289. else if (QR->size1 != b->size)
  290. {
  291. GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
  292. }
  293. else if (QR->size2 != x->size)
  294. {
  295. GSL_ERROR ("matrix size must match x size", GSL_EBADLEN);
  296. }
  297. else if (p->size != x->size)
  298. {
  299. GSL_ERROR ("permutation size must match x size", GSL_EBADLEN);
  300. }
  301. else
  302. {
  303. /* Copy x <- b */
  304. gsl_vector_memcpy (x, b);
  305. /* Solve R x = b, storing x inplace */
  306. gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, QR, x);
  307. gsl_permute_vector_inverse (p, x);
  308. return GSL_SUCCESS;
  309. }
  310. }
  311. int
  312. gsl_linalg_QRPT_Rsvx (const gsl_matrix * QR,
  313. const gsl_permutation * p,
  314. gsl_vector * x)
  315. {
  316. if (QR->size1 != QR->size2)
  317. {
  318. GSL_ERROR ("QR matrix must be square", GSL_ENOTSQR);
  319. }
  320. else if (QR->size2 != x->size)
  321. {
  322. GSL_ERROR ("matrix size must match x size", GSL_EBADLEN);
  323. }
  324. else if (p->size != x->size)
  325. {
  326. GSL_ERROR ("permutation size must match x size", GSL_EBADLEN);
  327. }
  328. else
  329. {
  330. /* Solve R x = b, storing x inplace */
  331. gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, QR, x);
  332. gsl_permute_vector_inverse (p, x);
  333. return GSL_SUCCESS;
  334. }
  335. }
  336. /* Update a Q R P^T factorisation for A P= Q R , A' = A + u v^T,
  337. Q' R' P^-1 = QR P^-1 + u v^T
  338. = Q (R + Q^T u v^T P ) P^-1
  339. = Q (R + w v^T P) P^-1
  340. where w = Q^T u.
  341. Algorithm from Golub and Van Loan, "Matrix Computations", Section
  342. 12.5 (Updating Matrix Factorizations, Rank-One Changes) */
  343. int
  344. gsl_linalg_QRPT_update (gsl_matrix * Q, gsl_matrix * R,
  345. const gsl_permutation * p,
  346. gsl_vector * w, const gsl_vector * v)
  347. {
  348. if (Q->size1 != Q->size2 || R->size1 != R->size2)
  349. {
  350. return GSL_ENOTSQR;
  351. }
  352. else if (R->size1 != Q->size2 || v->size != Q->size2 || w->size != Q->size2)
  353. {
  354. return GSL_EBADLEN;
  355. }
  356. else
  357. {
  358. size_t j, k;
  359. const size_t M = Q->size1;
  360. const size_t N = Q->size2;
  361. double w0;
  362. /* Apply Given's rotations to reduce w to (|w|, 0, 0, ... , 0)
  363. J_1^T .... J_(n-1)^T w = +/- |w| e_1
  364. simultaneously applied to R, H = J_1^T ... J^T_(n-1) R
  365. so that H is upper Hessenberg. (12.5.2) */
  366. for (k = N - 1; k > 0; k--)
  367. {
  368. double c, s;
  369. double wk = gsl_vector_get (w, k);
  370. double wkm1 = gsl_vector_get (w, k - 1);
  371. create_givens (wkm1, wk, &c, &s);
  372. apply_givens_vec (w, k - 1, k, c, s);
  373. apply_givens_qr (M, N, Q, R, k - 1, k, c, s);
  374. }
  375. w0 = gsl_vector_get (w, 0);
  376. /* Add in w v^T (Equation 12.5.3) */
  377. for (j = 0; j < N; j++)
  378. {
  379. double r0j = gsl_matrix_get (R, 0, j);
  380. size_t p_j = gsl_permutation_get (p, j);
  381. double vj = gsl_vector_get (v, p_j);
  382. gsl_matrix_set (R, 0, j, r0j + w0 * vj);
  383. }
  384. /* Apply Givens transformations R' = G_(n-1)^T ... G_1^T H
  385. Equation 12.5.4 */
  386. for (k = 1; k < N; k++)
  387. {
  388. double c, s;
  389. double diag = gsl_matrix_get (R, k - 1, k - 1);
  390. double offdiag = gsl_matrix_get (R, k, k - 1);
  391. create_givens (diag, offdiag, &c, &s);
  392. apply_givens_qr (M, N, Q, R, k - 1, k, c, s);
  393. }
  394. return GSL_SUCCESS;
  395. }
  396. }