gsl_linalg__lq.c 14 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568
  1. /* linalg/lq.c
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 Gerard Jungman, Brian Gough
  4. * Copyright (C) 2004 Joerg Wensch, modifications for LQ.
  5. *
  6. * This program is free software; you can redistribute it and/or modify
  7. * it under the terms of the GNU General Public License as published by
  8. * the Free Software Foundation; either version 3 of the License, or (at
  9. * your option) any later version.
  10. *
  11. * This program is distributed in the hope that it will be useful, but
  12. * WITHOUT ANY WARRANTY; without even the implied warranty of
  13. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
  14. * General Public License for more details.
  15. *
  16. * You should have received a copy of the GNU General Public License
  17. * along with this program; if not, write to the Free Software
  18. * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
  19. */
  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. /* Note: The standard in numerical linear algebra is to solve A x = b
  32. * resp. ||A x - b||_2 -> min by QR-decompositions where x, b are
  33. * column vectors.
  34. *
  35. * When the matrix A has a large number of rows it is much more
  36. * efficient to work with the transposed matrix A^T and to solve the
  37. * system x^T A = b^T resp. ||x^T A - b^T||_2 -> min. This is caused
  38. * by the row-oriented format in which GSL stores matrices. Therefore
  39. * the QR-decomposition of A has to be replaced by a LQ decomposition
  40. * of A^T
  41. *
  42. * The purpose of this package is to provide the algorithms to compute
  43. * the LQ-decomposition and to solve the linear equations resp. least
  44. * squares problems. The dimensions N, M of the matrix are switched
  45. * because here A will probably be a transposed matrix. We write x^T,
  46. * b^T,... for vectors the comments to emphasize that they are row
  47. * vectors.
  48. *
  49. * It may even be useful to transpose your matrix explicitly (assumed
  50. * that there are no memory restrictions) because this takes O(M x N)
  51. * computing time where the decompostion takes O(M x N^2) computing
  52. * time. */
  53. /* Factorise a general N x M matrix A into
  54. *
  55. * A = L Q
  56. *
  57. * where Q is orthogonal (M x M) and L is lower triangular (N x M).
  58. *
  59. * Q is stored as a packed set of Householder transformations in the
  60. * strict upper triangular part of the input matrix.
  61. *
  62. * R is stored in the diagonal and lower triangle of the input matrix.
  63. *
  64. * The full matrix for Q can be obtained as the product
  65. *
  66. * Q = Q_k .. Q_2 Q_1
  67. *
  68. * where k = MIN(M,N) and
  69. *
  70. * Q_i = (I - tau_i * v_i * v_i')
  71. *
  72. * and where v_i is a Householder vector
  73. *
  74. * v_i = [1, m(i+1,i), m(i+2,i), ... , m(M,i)]
  75. *
  76. * This storage scheme is the same as in LAPACK. */
  77. int
  78. gsl_linalg_LQ_decomp (gsl_matrix * A, gsl_vector * tau)
  79. {
  80. const size_t N = A->size1;
  81. const size_t M = A->size2;
  82. if (tau->size != GSL_MIN (M, N))
  83. {
  84. GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
  85. }
  86. else
  87. {
  88. size_t i;
  89. for (i = 0; i < GSL_MIN (M, N); i++)
  90. {
  91. /* Compute the Householder transformation to reduce the j-th
  92. column of the matrix to a multiple of the j-th unit vector */
  93. gsl_vector_view c_full = gsl_matrix_row (A, i);
  94. gsl_vector_view c = gsl_vector_subvector (&(c_full.vector), i, M-i);
  95. double tau_i = gsl_linalg_householder_transform (&(c.vector));
  96. gsl_vector_set (tau, i, tau_i);
  97. /* Apply the transformation to the remaining columns and
  98. update the norms */
  99. if (i + 1 < N)
  100. {
  101. gsl_matrix_view m = gsl_matrix_submatrix (A, i + 1, i, N - (i + 1), M - i );
  102. gsl_linalg_householder_mh (tau_i, &(c.vector), &(m.matrix));
  103. }
  104. }
  105. return GSL_SUCCESS;
  106. }
  107. }
  108. /* Solves the system x^T A = b^T using the LQ factorisation,
  109. * x^T L = b^T Q^T
  110. *
  111. * to obtain x. Based on SLATEC code.
  112. */
  113. int
  114. gsl_linalg_LQ_solve_T (const gsl_matrix * LQ, const gsl_vector * tau, const gsl_vector * b, gsl_vector * x)
  115. {
  116. if (LQ->size1 != LQ->size2)
  117. {
  118. GSL_ERROR ("LQ matrix must be square", GSL_ENOTSQR);
  119. }
  120. else if (LQ->size2 != b->size)
  121. {
  122. GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
  123. }
  124. else if (LQ->size1 != x->size)
  125. {
  126. GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
  127. }
  128. else
  129. {
  130. /* Copy x <- b */
  131. gsl_vector_memcpy (x, b);
  132. /* Solve for x */
  133. gsl_linalg_LQ_svx_T (LQ, tau, x);
  134. return GSL_SUCCESS;
  135. }
  136. }
  137. /* Solves the system x^T A = b^T in place using the LQ factorisation,
  138. *
  139. * x^T L = b^T Q^T
  140. *
  141. * to obtain x. Based on SLATEC code.
  142. */
  143. int
  144. gsl_linalg_LQ_svx_T (const gsl_matrix * LQ, const gsl_vector * tau, gsl_vector * x)
  145. {
  146. if (LQ->size1 != LQ->size2)
  147. {
  148. GSL_ERROR ("LQ matrix must be square", GSL_ENOTSQR);
  149. }
  150. else if (LQ->size1 != x->size)
  151. {
  152. GSL_ERROR ("matrix size must match x/rhs size", GSL_EBADLEN);
  153. }
  154. else
  155. {
  156. /* compute rhs = Q^T b */
  157. gsl_linalg_LQ_vecQT (LQ, tau, x);
  158. /* Solve R x = rhs, storing x in-place */
  159. gsl_blas_dtrsv (CblasLower, CblasTrans, CblasNonUnit, LQ, x);
  160. return GSL_SUCCESS;
  161. }
  162. }
  163. /* Find the least squares solution to the overdetermined system
  164. *
  165. * x^T A = b^T
  166. *
  167. * for M >= N using the LQ factorization A = L Q.
  168. */
  169. int
  170. gsl_linalg_LQ_lssolve_T (const gsl_matrix * LQ, const gsl_vector * tau, const gsl_vector * b, gsl_vector * x, gsl_vector * residual)
  171. {
  172. const size_t N = LQ->size1;
  173. const size_t M = LQ->size2;
  174. if (M < N)
  175. {
  176. GSL_ERROR ("LQ matrix must have M>=N", GSL_EBADLEN);
  177. }
  178. else if (M != b->size)
  179. {
  180. GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
  181. }
  182. else if (N != x->size)
  183. {
  184. GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
  185. }
  186. else if (M != residual->size)
  187. {
  188. GSL_ERROR ("matrix size must match residual size", GSL_EBADLEN);
  189. }
  190. else
  191. {
  192. gsl_matrix_const_view L = gsl_matrix_const_submatrix (LQ, 0, 0, N, N);
  193. gsl_vector_view c = gsl_vector_subvector(residual, 0, N);
  194. gsl_vector_memcpy(residual, b);
  195. /* compute rhs = b^T Q^T */
  196. gsl_linalg_LQ_vecQT (LQ, tau, residual);
  197. /* Solve x^T L = rhs */
  198. gsl_vector_memcpy(x, &(c.vector));
  199. gsl_blas_dtrsv (CblasLower, CblasTrans, CblasNonUnit, &(L.matrix), x);
  200. /* Compute residual = b^T - x^T A = (b^T Q^T - x^T L) Q */
  201. gsl_vector_set_zero(&(c.vector));
  202. gsl_linalg_LQ_vecQ(LQ, tau, residual);
  203. return GSL_SUCCESS;
  204. }
  205. }
  206. int
  207. gsl_linalg_LQ_Lsolve_T (const gsl_matrix * LQ, const gsl_vector * b, gsl_vector * x)
  208. {
  209. if (LQ->size1 != LQ->size2)
  210. {
  211. GSL_ERROR ("LQ matrix must be square", GSL_ENOTSQR);
  212. }
  213. else if (LQ->size1 != b->size)
  214. {
  215. GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
  216. }
  217. else if (LQ->size1 != x->size)
  218. {
  219. GSL_ERROR ("matrix size must match x size", GSL_EBADLEN);
  220. }
  221. else
  222. {
  223. /* Copy x <- b */
  224. gsl_vector_memcpy (x, b);
  225. /* Solve R x = b, storing x in-place */
  226. gsl_blas_dtrsv (CblasLower, CblasTrans, CblasNonUnit, LQ, x);
  227. return GSL_SUCCESS;
  228. }
  229. }
  230. int
  231. gsl_linalg_LQ_Lsvx_T (const gsl_matrix * LQ, gsl_vector * x)
  232. {
  233. if (LQ->size1 != LQ->size2)
  234. {
  235. GSL_ERROR ("LQ matrix must be square", GSL_ENOTSQR);
  236. }
  237. else if (LQ->size2 != x->size)
  238. {
  239. GSL_ERROR ("matrix size must match rhs size", GSL_EBADLEN);
  240. }
  241. else
  242. {
  243. /* Solve x^T L = b^T, storing x in-place */
  244. gsl_blas_dtrsv (CblasLower, CblasTrans, CblasNonUnit, LQ, x);
  245. return GSL_SUCCESS;
  246. }
  247. }
  248. int
  249. gsl_linalg_L_solve_T (const gsl_matrix * L, const gsl_vector * b, gsl_vector * x)
  250. {
  251. if (L->size1 != L->size2)
  252. {
  253. GSL_ERROR ("R matrix must be square", GSL_ENOTSQR);
  254. }
  255. else if (L->size2 != b->size)
  256. {
  257. GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
  258. }
  259. else if (L->size1 != x->size)
  260. {
  261. GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
  262. }
  263. else
  264. {
  265. /* Copy x <- b */
  266. gsl_vector_memcpy (x, b);
  267. /* Solve R x = b, storing x inplace in b */
  268. gsl_blas_dtrsv (CblasLower, CblasTrans, CblasNonUnit, L, x);
  269. return GSL_SUCCESS;
  270. }
  271. }
  272. int
  273. gsl_linalg_LQ_vecQT (const gsl_matrix * LQ, const gsl_vector * tau, gsl_vector * v)
  274. {
  275. const size_t N = LQ->size1;
  276. const size_t M = LQ->size2;
  277. if (tau->size != GSL_MIN (M, N))
  278. {
  279. GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
  280. }
  281. else if (v->size != M)
  282. {
  283. GSL_ERROR ("vector size must be M", GSL_EBADLEN);
  284. }
  285. else
  286. {
  287. size_t i;
  288. /* compute v Q^T */
  289. for (i = 0; i < GSL_MIN (M, N); i++)
  290. {
  291. gsl_vector_const_view c = gsl_matrix_const_row (LQ, i);
  292. gsl_vector_const_view h = gsl_vector_const_subvector (&(c.vector),
  293. i, M - i);
  294. gsl_vector_view w = gsl_vector_subvector (v, i, M - i);
  295. double ti = gsl_vector_get (tau, i);
  296. gsl_linalg_householder_hv (ti, &(h.vector), &(w.vector));
  297. }
  298. return GSL_SUCCESS;
  299. }
  300. }
  301. int
  302. gsl_linalg_LQ_vecQ (const gsl_matrix * LQ, const gsl_vector * tau, gsl_vector * v)
  303. {
  304. const size_t N = LQ->size1;
  305. const size_t M = LQ->size2;
  306. if (tau->size != GSL_MIN (M, N))
  307. {
  308. GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
  309. }
  310. else if (v->size != M)
  311. {
  312. GSL_ERROR ("vector size must be M", GSL_EBADLEN);
  313. }
  314. else
  315. {
  316. size_t i;
  317. /* compute v Q^T */
  318. for (i = GSL_MIN (M, N); i > 0 && i--;)
  319. {
  320. gsl_vector_const_view c = gsl_matrix_const_row (LQ, i);
  321. gsl_vector_const_view h = gsl_vector_const_subvector (&(c.vector),
  322. i, M - i);
  323. gsl_vector_view w = gsl_vector_subvector (v, i, M - i);
  324. double ti = gsl_vector_get (tau, i);
  325. gsl_linalg_householder_hv (ti, &(h.vector), &(w.vector));
  326. }
  327. return GSL_SUCCESS;
  328. }
  329. }
  330. /* Form the orthogonal matrix Q from the packed LQ matrix */
  331. int
  332. gsl_linalg_LQ_unpack (const gsl_matrix * LQ, const gsl_vector * tau, gsl_matrix * Q, gsl_matrix * L)
  333. {
  334. const size_t N = LQ->size1;
  335. const size_t M = LQ->size2;
  336. if (Q->size1 != M || Q->size2 != M)
  337. {
  338. GSL_ERROR ("Q matrix must be M x M", GSL_ENOTSQR);
  339. }
  340. else if (L->size1 != N || L->size2 != M)
  341. {
  342. GSL_ERROR ("R matrix must be N x M", GSL_ENOTSQR);
  343. }
  344. else if (tau->size != GSL_MIN (M, N))
  345. {
  346. GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
  347. }
  348. else
  349. {
  350. size_t i, j, l_border;
  351. /* Initialize Q to the identity */
  352. gsl_matrix_set_identity (Q);
  353. for (i = GSL_MIN (M, N); i > 0 && i--;)
  354. {
  355. gsl_vector_const_view c = gsl_matrix_const_row (LQ, i);
  356. gsl_vector_const_view h = gsl_vector_const_subvector (&c.vector,
  357. i, M - i);
  358. gsl_matrix_view m = gsl_matrix_submatrix (Q, i, i, M - i, M - i);
  359. double ti = gsl_vector_get (tau, i);
  360. gsl_linalg_householder_mh (ti, &h.vector, &m.matrix);
  361. }
  362. /* Form the lower triangular matrix L from a packed LQ matrix */
  363. for (i = 0; i < N; i++)
  364. {
  365. l_border=GSL_MIN(i,M-1);
  366. for (j = 0; j <= l_border ; j++)
  367. gsl_matrix_set (L, i, j, gsl_matrix_get (LQ, i, j));
  368. for (j = l_border+1; j < M; j++)
  369. gsl_matrix_set (L, i, j, 0.0);
  370. }
  371. return GSL_SUCCESS;
  372. }
  373. }
  374. /* Update a LQ factorisation for A= L Q , A' = A + v u^T,
  375. * L' Q' = LQ + v u^T
  376. * = (L + v u^T Q^T) Q
  377. * = (L + v w^T) Q
  378. *
  379. * where w = Q u.
  380. *
  381. * Algorithm from Golub and Van Loan, "Matrix Computations", Section
  382. * 12.5 (Updating Matrix Factorizations, Rank-One Changes)
  383. */
  384. int
  385. gsl_linalg_LQ_update (gsl_matrix * Q, gsl_matrix * L,
  386. const gsl_vector * v, gsl_vector * w)
  387. {
  388. const size_t N = L->size1;
  389. const size_t M = L->size2;
  390. if (Q->size1 != M || Q->size2 != M)
  391. {
  392. GSL_ERROR ("Q matrix must be N x N if L is M x N", GSL_ENOTSQR);
  393. }
  394. else if (w->size != M)
  395. {
  396. GSL_ERROR ("w must be length N if L is M x N", GSL_EBADLEN);
  397. }
  398. else if (v->size != N)
  399. {
  400. GSL_ERROR ("v must be length M if L is M x N", GSL_EBADLEN);
  401. }
  402. else
  403. {
  404. size_t j, k;
  405. double w0;
  406. /* Apply Given's rotations to reduce w to (|w|, 0, 0, ... , 0)
  407. J_1^T .... J_(n-1)^T w = +/- |w| e_1
  408. simultaneously applied to L, H = J_1^T ... J^T_(n-1) L
  409. so that H is upper Hessenberg. (12.5.2) */
  410. for (k = M - 1; k > 0; k--)
  411. {
  412. double c, s;
  413. double wk = gsl_vector_get (w, k);
  414. double wkm1 = gsl_vector_get (w, k - 1);
  415. create_givens (wkm1, wk, &c, &s);
  416. apply_givens_vec (w, k - 1, k, c, s);
  417. apply_givens_lq (M, N, Q, L, k - 1, k, c, s);
  418. }
  419. w0 = gsl_vector_get (w, 0);
  420. /* Add in v w^T (Equation 12.5.3) */
  421. for (j = 0; j < N; j++)
  422. {
  423. double lj0 = gsl_matrix_get (L, j, 0);
  424. double vj = gsl_vector_get (v, j);
  425. gsl_matrix_set (L, j, 0, lj0 + w0 * vj);
  426. }
  427. /* Apply Givens transformations L' = G_(n-1)^T ... G_1^T H
  428. Equation 12.5.4 */
  429. for (k = 1; k < GSL_MIN(M,N+1); k++)
  430. {
  431. double c, s;
  432. double diag = gsl_matrix_get (L, k - 1, k - 1);
  433. double offdiag = gsl_matrix_get (L, k - 1 , k);
  434. create_givens (diag, offdiag, &c, &s);
  435. apply_givens_lq (M, N, Q, L, k - 1, k, c, s);
  436. gsl_matrix_set (L, k - 1, k, 0.0); /* exact zero of G^T */
  437. }
  438. return GSL_SUCCESS;
  439. }
  440. }
  441. int
  442. gsl_linalg_LQ_LQsolve (gsl_matrix * Q, gsl_matrix * L, const gsl_vector * b, gsl_vector * x)
  443. {
  444. const size_t N = L->size1;
  445. const size_t M = L->size2;
  446. if (M != N)
  447. {
  448. return GSL_ENOTSQR;
  449. }
  450. else if (Q->size1 != M || b->size != M || x->size != M)
  451. {
  452. return GSL_EBADLEN;
  453. }
  454. else
  455. {
  456. /* compute sol = b^T Q^T */
  457. gsl_blas_dgemv (CblasNoTrans, 1.0, Q, b, 0.0, x);
  458. /* Solve x^T L = sol, storing x in-place */
  459. gsl_blas_dtrsv (CblasLower, CblasTrans, CblasNonUnit, L, x);
  460. return GSL_SUCCESS;
  461. }
  462. }