gsl_linalg__bidiag.c 9.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365
  1. /* linalg/bidiag.c
  2. *
  3. * Copyright (C) 2001, 2007 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. /* Factorise a matrix A into
  20. *
  21. * A = U B V^T
  22. *
  23. * where U and V are orthogonal and B is upper bidiagonal.
  24. *
  25. * On exit, B is stored in the diagonal and first superdiagonal of A.
  26. *
  27. * U is stored as a packed set of Householder transformations in the
  28. * lower triangular part of the input matrix below the diagonal.
  29. *
  30. * V is stored as a packed set of Householder transformations in the
  31. * upper triangular part of the input matrix above the first
  32. * superdiagonal.
  33. *
  34. * The full matrix for U can be obtained as the product
  35. *
  36. * U = U_1 U_2 .. U_N
  37. *
  38. * where
  39. *
  40. * U_i = (I - tau_i * u_i * u_i')
  41. *
  42. * and where u_i is a Householder vector
  43. *
  44. * u_i = [0, .. , 0, 1, A(i+1,i), A(i+3,i), .. , A(M,i)]
  45. *
  46. * The full matrix for V can be obtained as the product
  47. *
  48. * V = V_1 V_2 .. V_(N-2)
  49. *
  50. * where
  51. *
  52. * V_i = (I - tau_i * v_i * v_i')
  53. *
  54. * and where v_i is a Householder vector
  55. *
  56. * v_i = [0, .. , 0, 1, A(i,i+2), A(i,i+3), .. , A(i,N)]
  57. *
  58. * See Golub & Van Loan, "Matrix Computations" (3rd ed), Algorithm 5.4.2
  59. *
  60. * Note: this description uses 1-based indices. The code below uses
  61. * 0-based indices
  62. */
  63. #include "gsl__config.h"
  64. #include <stdlib.h>
  65. #include "gsl_math.h"
  66. #include "gsl_vector.h"
  67. #include "gsl_matrix.h"
  68. #include "gsl_blas.h"
  69. #include "gsl_linalg.h"
  70. int
  71. gsl_linalg_bidiag_decomp (gsl_matrix * A, gsl_vector * tau_U, gsl_vector * tau_V)
  72. {
  73. if (A->size1 < A->size2)
  74. {
  75. GSL_ERROR ("bidiagonal decomposition requires M>=N", GSL_EBADLEN);
  76. }
  77. else if (tau_U->size != A->size2)
  78. {
  79. GSL_ERROR ("size of tau_U must be N", GSL_EBADLEN);
  80. }
  81. else if (tau_V->size + 1 != A->size2)
  82. {
  83. GSL_ERROR ("size of tau_V must be (N - 1)", GSL_EBADLEN);
  84. }
  85. else
  86. {
  87. const size_t M = A->size1;
  88. const size_t N = A->size2;
  89. size_t i;
  90. for (i = 0 ; i < N; i++)
  91. {
  92. /* Apply Householder transformation to current column */
  93. {
  94. gsl_vector_view c = gsl_matrix_column (A, i);
  95. gsl_vector_view v = gsl_vector_subvector (&c.vector, i, M - i);
  96. double tau_i = gsl_linalg_householder_transform (&v.vector);
  97. /* Apply the transformation to the remaining columns */
  98. if (i + 1 < N)
  99. {
  100. gsl_matrix_view m =
  101. gsl_matrix_submatrix (A, i, i + 1, M - i, N - (i + 1));
  102. gsl_linalg_householder_hm (tau_i, &v.vector, &m.matrix);
  103. }
  104. gsl_vector_set (tau_U, i, tau_i);
  105. }
  106. /* Apply Householder transformation to current row */
  107. if (i + 1 < N)
  108. {
  109. gsl_vector_view r = gsl_matrix_row (A, i);
  110. gsl_vector_view v = gsl_vector_subvector (&r.vector, i + 1, N - (i + 1));
  111. double tau_i = gsl_linalg_householder_transform (&v.vector);
  112. /* Apply the transformation to the remaining rows */
  113. if (i + 1 < M)
  114. {
  115. gsl_matrix_view m =
  116. gsl_matrix_submatrix (A, i+1, i+1, M - (i+1), N - (i+1));
  117. gsl_linalg_householder_mh (tau_i, &v.vector, &m.matrix);
  118. }
  119. gsl_vector_set (tau_V, i, tau_i);
  120. }
  121. }
  122. }
  123. return GSL_SUCCESS;
  124. }
  125. /* Form the orthogonal matrices U, V, diagonal d and superdiagonal sd
  126. from the packed bidiagonal matrix A */
  127. int
  128. gsl_linalg_bidiag_unpack (const gsl_matrix * A,
  129. const gsl_vector * tau_U,
  130. gsl_matrix * U,
  131. const gsl_vector * tau_V,
  132. gsl_matrix * V,
  133. gsl_vector * diag,
  134. gsl_vector * superdiag)
  135. {
  136. const size_t M = A->size1;
  137. const size_t N = A->size2;
  138. const size_t K = GSL_MIN(M, N);
  139. if (M < N)
  140. {
  141. GSL_ERROR ("matrix A must have M >= N", GSL_EBADLEN);
  142. }
  143. else if (tau_U->size != K)
  144. {
  145. GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
  146. }
  147. else if (tau_V->size + 1 != K)
  148. {
  149. GSL_ERROR ("size of tau must be MIN(M,N) - 1", GSL_EBADLEN);
  150. }
  151. else if (U->size1 != M || U->size2 != N)
  152. {
  153. GSL_ERROR ("size of U must be M x N", GSL_EBADLEN);
  154. }
  155. else if (V->size1 != N || V->size2 != N)
  156. {
  157. GSL_ERROR ("size of V must be N x N", GSL_EBADLEN);
  158. }
  159. else if (diag->size != K)
  160. {
  161. GSL_ERROR ("size of diagonal must match size of A", GSL_EBADLEN);
  162. }
  163. else if (superdiag->size + 1 != K)
  164. {
  165. GSL_ERROR ("size of subdiagonal must be (diagonal size - 1)", GSL_EBADLEN);
  166. }
  167. else
  168. {
  169. size_t i, j;
  170. /* Copy diagonal into diag */
  171. for (i = 0; i < N; i++)
  172. {
  173. double Aii = gsl_matrix_get (A, i, i);
  174. gsl_vector_set (diag, i, Aii);
  175. }
  176. /* Copy superdiagonal into superdiag */
  177. for (i = 0; i < N - 1; i++)
  178. {
  179. double Aij = gsl_matrix_get (A, i, i+1);
  180. gsl_vector_set (superdiag, i, Aij);
  181. }
  182. /* Initialize V to the identity */
  183. gsl_matrix_set_identity (V);
  184. for (i = N - 1; i > 0 && i--;)
  185. {
  186. /* Householder row transformation to accumulate V */
  187. gsl_vector_const_view r = gsl_matrix_const_row (A, i);
  188. gsl_vector_const_view h =
  189. gsl_vector_const_subvector (&r.vector, i + 1, N - (i+1));
  190. double ti = gsl_vector_get (tau_V, i);
  191. gsl_matrix_view m =
  192. gsl_matrix_submatrix (V, i + 1, i + 1, N-(i+1), N-(i+1));
  193. gsl_linalg_householder_hm (ti, &h.vector, &m.matrix);
  194. }
  195. /* Initialize U to the identity */
  196. gsl_matrix_set_identity (U);
  197. for (j = N; j > 0 && j--;)
  198. {
  199. /* Householder column transformation to accumulate U */
  200. gsl_vector_const_view c = gsl_matrix_const_column (A, j);
  201. gsl_vector_const_view h = gsl_vector_const_subvector (&c.vector, j, M - j);
  202. double tj = gsl_vector_get (tau_U, j);
  203. gsl_matrix_view m =
  204. gsl_matrix_submatrix (U, j, j, M-j, N-j);
  205. gsl_linalg_householder_hm (tj, &h.vector, &m.matrix);
  206. }
  207. return GSL_SUCCESS;
  208. }
  209. }
  210. int
  211. gsl_linalg_bidiag_unpack2 (gsl_matrix * A,
  212. gsl_vector * tau_U,
  213. gsl_vector * tau_V,
  214. gsl_matrix * V)
  215. {
  216. const size_t M = A->size1;
  217. const size_t N = A->size2;
  218. const size_t K = GSL_MIN(M, N);
  219. if (M < N)
  220. {
  221. GSL_ERROR ("matrix A must have M >= N", GSL_EBADLEN);
  222. }
  223. else if (tau_U->size != K)
  224. {
  225. GSL_ERROR ("size of tau must be MIN(M,N)", GSL_EBADLEN);
  226. }
  227. else if (tau_V->size + 1 != K)
  228. {
  229. GSL_ERROR ("size of tau must be MIN(M,N) - 1", GSL_EBADLEN);
  230. }
  231. else if (V->size1 != N || V->size2 != N)
  232. {
  233. GSL_ERROR ("size of V must be N x N", GSL_EBADLEN);
  234. }
  235. else
  236. {
  237. size_t i, j;
  238. /* Initialize V to the identity */
  239. gsl_matrix_set_identity (V);
  240. for (i = N - 1; i > 0 && i--;)
  241. {
  242. /* Householder row transformation to accumulate V */
  243. gsl_vector_const_view r = gsl_matrix_const_row (A, i);
  244. gsl_vector_const_view h =
  245. gsl_vector_const_subvector (&r.vector, i + 1, N - (i+1));
  246. double ti = gsl_vector_get (tau_V, i);
  247. gsl_matrix_view m =
  248. gsl_matrix_submatrix (V, i + 1, i + 1, N-(i+1), N-(i+1));
  249. gsl_linalg_householder_hm (ti, &h.vector, &m.matrix);
  250. }
  251. /* Copy superdiagonal into tau_v */
  252. for (i = 0; i < N - 1; i++)
  253. {
  254. double Aij = gsl_matrix_get (A, i, i+1);
  255. gsl_vector_set (tau_V, i, Aij);
  256. }
  257. /* Allow U to be unpacked into the same memory as A, copy
  258. diagonal into tau_U */
  259. for (j = N; j > 0 && j--;)
  260. {
  261. /* Householder column transformation to accumulate U */
  262. double tj = gsl_vector_get (tau_U, j);
  263. double Ajj = gsl_matrix_get (A, j, j);
  264. gsl_matrix_view m = gsl_matrix_submatrix (A, j, j, M-j, N-j);
  265. gsl_vector_set (tau_U, j, Ajj);
  266. gsl_linalg_householder_hm1 (tj, &m.matrix);
  267. }
  268. return GSL_SUCCESS;
  269. }
  270. }
  271. int
  272. gsl_linalg_bidiag_unpack_B (const gsl_matrix * A,
  273. gsl_vector * diag,
  274. gsl_vector * superdiag)
  275. {
  276. const size_t M = A->size1;
  277. const size_t N = A->size2;
  278. const size_t K = GSL_MIN(M, N);
  279. if (diag->size != K)
  280. {
  281. GSL_ERROR ("size of diagonal must match size of A", GSL_EBADLEN);
  282. }
  283. else if (superdiag->size + 1 != K)
  284. {
  285. GSL_ERROR ("size of subdiagonal must be (matrix size - 1)", GSL_EBADLEN);
  286. }
  287. else
  288. {
  289. size_t i;
  290. /* Copy diagonal into diag */
  291. for (i = 0; i < K; i++)
  292. {
  293. double Aii = gsl_matrix_get (A, i, i);
  294. gsl_vector_set (diag, i, Aii);
  295. }
  296. /* Copy superdiagonal into superdiag */
  297. for (i = 0; i < K - 1; i++)
  298. {
  299. double Aij = gsl_matrix_get (A, i, i+1);
  300. gsl_vector_set (superdiag, i, Aij);
  301. }
  302. return GSL_SUCCESS;
  303. }
  304. }