gsl_linalg__hessenberg.c 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434
  1. /* linalg/hessenberg.c
  2. *
  3. * Copyright (C) 2006 Patrick Alken
  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 "gsl_linalg.h"
  21. #include "gsl_matrix.h"
  22. #include "gsl_vector.h"
  23. /*
  24. gsl_linalg_hessenberg_decomp()
  25. Compute the Householder reduction to Hessenberg form of a
  26. square N-by-N matrix A.
  27. H = U^t A U
  28. See Golub & Van Loan, "Matrix Computations" (3rd ed), algorithm
  29. 7.4.2
  30. Inputs: A - matrix to reduce
  31. tau - where to store scalar factors in Householder
  32. matrices; this vector must be of length N,
  33. where N is the order of A
  34. Return: GSL_SUCCESS unless error occurs
  35. Notes: on output, the upper triangular portion of A (including
  36. the diagaonal and subdiagonal) contains the Hessenberg matrix.
  37. The lower triangular portion (below the subdiagonal) contains
  38. the Householder vectors which can be used to construct
  39. the similarity transform matrix U.
  40. The matrix U is
  41. U = U(1) U(2) ... U(n - 2)
  42. where
  43. U(i) = I - tau(i) * v(i) * v(i)^t
  44. and the vector v(i) is stored in column i of the matrix A
  45. underneath the subdiagonal. So the first element of v(i)
  46. is stored in row i + 2, column i, the second element at
  47. row i + 3, column i, and so on.
  48. Also note that for the purposes of computing U(i),
  49. v(1:i) = 0, v(i + 1) = 1, and v(i+2:n) is what is stored in
  50. column i of A beneath the subdiagonal.
  51. */
  52. int
  53. gsl_linalg_hessenberg_decomp(gsl_matrix *A, gsl_vector *tau)
  54. {
  55. const size_t N = A->size1;
  56. if (N != A->size2)
  57. {
  58. GSL_ERROR ("Hessenberg reduction requires square matrix",
  59. GSL_ENOTSQR);
  60. }
  61. else if (N != tau->size)
  62. {
  63. GSL_ERROR ("tau vector must match matrix size", GSL_EBADLEN);
  64. }
  65. else if (N < 3)
  66. {
  67. /* nothing to do */
  68. return GSL_SUCCESS;
  69. }
  70. else
  71. {
  72. size_t i; /* looping */
  73. gsl_vector_view c, /* matrix column */
  74. hv; /* householder vector */
  75. gsl_matrix_view m;
  76. double tau_i; /* beta in algorithm 7.4.2 */
  77. for (i = 0; i < N - 2; ++i)
  78. {
  79. /*
  80. * make a copy of A(i + 1:n, i) and store it in the section
  81. * of 'tau' that we haven't stored coefficients in yet
  82. */
  83. c = gsl_matrix_subcolumn(A, i, i + 1, N - i - 1);
  84. hv = gsl_vector_subvector(tau, i + 1, N - (i + 1));
  85. gsl_vector_memcpy(&hv.vector, &c.vector);
  86. /* compute householder transformation of A(i+1:n,i) */
  87. tau_i = gsl_linalg_householder_transform(&hv.vector);
  88. /* apply left householder matrix (I - tau_i v v') to A */
  89. m = gsl_matrix_submatrix(A, i + 1, i, N - (i + 1), N - i);
  90. gsl_linalg_householder_hm(tau_i, &hv.vector, &m.matrix);
  91. /* apply right householder matrix (I - tau_i v v') to A */
  92. m = gsl_matrix_submatrix(A, 0, i + 1, N, N - (i + 1));
  93. gsl_linalg_householder_mh(tau_i, &hv.vector, &m.matrix);
  94. /* save Householder coefficient */
  95. gsl_vector_set(tau, i, tau_i);
  96. /*
  97. * store Householder vector below the subdiagonal in column
  98. * i of the matrix. hv(1) does not need to be stored since
  99. * it is always 1.
  100. */
  101. c = gsl_vector_subvector(&c.vector, 1, c.vector.size - 1);
  102. hv = gsl_vector_subvector(&hv.vector, 1, hv.vector.size - 1);
  103. gsl_vector_memcpy(&c.vector, &hv.vector);
  104. }
  105. return GSL_SUCCESS;
  106. }
  107. } /* gsl_linalg_hessenberg_decomp() */
  108. /*
  109. gsl_linalg_hessenberg_unpack()
  110. Construct the matrix U which transforms a matrix A into
  111. its upper Hessenberg form:
  112. H = U^t A U
  113. by unpacking the information stored in H from gsl_linalg_hessenberg().
  114. U is a product of Householder matrices:
  115. U = U(1) U(2) ... U(n - 2)
  116. where
  117. U(i) = I - tau(i) * v(i) * v(i)^t
  118. The v(i) are stored in the lower triangular part of H by
  119. gsl_linalg_hessenberg(). The tau(i) are stored in the vector tau.
  120. Inputs: H - Hessenberg matrix computed from
  121. gsl_linalg_hessenberg()
  122. tau - tau vector computed from gsl_linalg_hessenberg()
  123. U - (output) where to store similarity matrix
  124. Return: success or error
  125. */
  126. int
  127. gsl_linalg_hessenberg_unpack(gsl_matrix * H, gsl_vector * tau,
  128. gsl_matrix * U)
  129. {
  130. int s;
  131. gsl_matrix_set_identity(U);
  132. s = gsl_linalg_hessenberg_unpack_accum(H, tau, U);
  133. return s;
  134. } /* gsl_linalg_hessenberg_unpack() */
  135. /*
  136. gsl_linalg_hessenberg_unpack_accum()
  137. This routine is the same as gsl_linalg_hessenberg_unpack(), except
  138. instead of storing the similarity matrix in U, it accumulates it,
  139. so that
  140. U -> U * [ U(1) U(2) ... U(n - 2) ]
  141. instead of:
  142. U -> U(1) U(2) ... U(n - 2)
  143. Inputs: H - Hessenberg matrix computed from
  144. gsl_linalg_hessenberg()
  145. tau - tau vector computed from gsl_linalg_hessenberg()
  146. V - (input/output) where to accumulate similarity matrix
  147. Return: success or error
  148. Notes: 1) On input, V needs to be initialized. The Householder matrices
  149. are accumulated into V, so on output,
  150. V_out = V_in * U(1) * U(2) * ... * U(n - 2)
  151. so if you just want the product of the Householder matrices,
  152. initialize V to the identity matrix before calling this
  153. function.
  154. 2) V does not have to be square, but must have the same
  155. number of columns as the order of H
  156. */
  157. int
  158. gsl_linalg_hessenberg_unpack_accum(gsl_matrix * H, gsl_vector * tau,
  159. gsl_matrix * V)
  160. {
  161. const size_t N = H->size1;
  162. if (N != H->size2)
  163. {
  164. GSL_ERROR ("Hessenberg reduction requires square matrix",
  165. GSL_ENOTSQR);
  166. }
  167. else if (N != tau->size)
  168. {
  169. GSL_ERROR ("tau vector must match matrix size", GSL_EBADLEN);
  170. }
  171. else if (N != V->size2)
  172. {
  173. GSL_ERROR ("V matrix has wrong dimension", GSL_EBADLEN);
  174. }
  175. else
  176. {
  177. size_t j; /* looping */
  178. double tau_j; /* householder coefficient */
  179. gsl_vector_view c, /* matrix column */
  180. hv; /* householder vector */
  181. gsl_matrix_view m;
  182. if (N < 3)
  183. {
  184. /* nothing to do */
  185. return GSL_SUCCESS;
  186. }
  187. for (j = 0; j < (N - 2); ++j)
  188. {
  189. c = gsl_matrix_column(H, j);
  190. tau_j = gsl_vector_get(tau, j);
  191. /*
  192. * get a view to the householder vector in column j, but
  193. * make sure hv(2) starts at the element below the
  194. * subdiagonal, since hv(1) was never stored and is always
  195. * 1
  196. */
  197. hv = gsl_vector_subvector(&c.vector, j + 1, N - (j + 1));
  198. /*
  199. * Only operate on part of the matrix since the first
  200. * j + 1 entries of the real householder vector are 0
  201. *
  202. * V -> V * U(j)
  203. *
  204. * Note here that V->size1 is not necessarily equal to N
  205. */
  206. m = gsl_matrix_submatrix(V, 0, j + 1, V->size1, N - (j + 1));
  207. /* apply right Householder matrix to V */
  208. gsl_linalg_householder_mh(tau_j, &hv.vector, &m.matrix);
  209. }
  210. return GSL_SUCCESS;
  211. }
  212. } /* gsl_linalg_hessenberg_unpack_accum() */
  213. /*
  214. gsl_linalg_hessenberg_set_zero()
  215. Zero out the lower triangular portion of the Hessenberg matrix H.
  216. This is useful when Householder vectors may be stored in the lower
  217. part of H, but eigenvalue solvers need some scratch space with zeros.
  218. */
  219. int
  220. gsl_linalg_hessenberg_set_zero(gsl_matrix * H)
  221. {
  222. const size_t N = H->size1;
  223. size_t i, j;
  224. if (N < 3)
  225. return GSL_SUCCESS;
  226. for (j = 0; j < N - 2; ++j)
  227. {
  228. for (i = j + 2; i < N; ++i)
  229. {
  230. gsl_matrix_set(H, i, j, 0.0);
  231. }
  232. }
  233. return GSL_SUCCESS;
  234. } /* gsl_linalg_hessenberg_set_zero() */
  235. /*
  236. gsl_linalg_hessenberg_submatrix()
  237. This routine does the same thing as gsl_linalg_hessenberg(),
  238. except that it operates on a submatrix of a larger matrix, but
  239. updates the larger matrix with the Householder transformations.
  240. For example, suppose
  241. M = [ M_{11} | M_{12} | M_{13} ]
  242. [ 0 | A | M_{23} ]
  243. [ 0 | 0 | M_{33} ]
  244. where M_{11} and M_{33} are already in Hessenberg form, and we
  245. just want to reduce A to Hessenberg form. Applying the transformations
  246. to A alone will cause the larger matrix M to lose its similarity
  247. information. So this routine updates M_{12} and M_{23} as A gets
  248. reduced.
  249. Inputs: M - total matrix
  250. A - (sub)matrix to reduce
  251. top - row index of top of A in M
  252. tau - where to store scalar factors in Householder
  253. matrices; this vector must be of length N,
  254. where N is the order of A
  255. Return: GSL_SUCCESS unless error occurs
  256. Notes: on output, the upper triangular portion of A (including
  257. the diagaonal and subdiagonal) contains the Hessenberg matrix.
  258. The lower triangular portion (below the subdiagonal) contains
  259. the Householder vectors which can be used to construct
  260. the similarity transform matrix U.
  261. The matrix U is
  262. U = U(1) U(2) ... U(n - 2)
  263. where
  264. U(i) = I - tau(i) * v(i) * v(i)^t
  265. and the vector v(i) is stored in column i of the matrix A
  266. underneath the subdiagonal. So the first element of v(i)
  267. is stored in row i + 2, column i, the second element at
  268. row i + 3, column i, and so on.
  269. Also note that for the purposes of computing U(i),
  270. v(1:i) = 0, v(i + 1) = 1, and v(i+2:n) is what is stored in
  271. column i of A beneath the subdiagonal.
  272. */
  273. int
  274. gsl_linalg_hessenberg_submatrix(gsl_matrix *M, gsl_matrix *A, size_t top,
  275. gsl_vector *tau)
  276. {
  277. const size_t N = A->size1;
  278. const size_t N_M = M->size1;
  279. if (N != A->size2)
  280. {
  281. GSL_ERROR ("Hessenberg reduction requires square matrix",
  282. GSL_ENOTSQR);
  283. }
  284. else if (N != tau->size)
  285. {
  286. GSL_ERROR ("tau vector must match matrix size", GSL_EBADLEN);
  287. }
  288. else if (N < 3)
  289. {
  290. /* nothing to do */
  291. return GSL_SUCCESS;
  292. }
  293. else
  294. {
  295. size_t i; /* looping */
  296. gsl_vector_view c, /* matrix column */
  297. hv; /* householder vector */
  298. gsl_matrix_view m;
  299. double tau_i; /* beta in algorithm 7.4.2 */
  300. for (i = 0; i < N - 2; ++i)
  301. {
  302. /*
  303. * make a copy of A(i + 1:n, i) and store it in the section
  304. * of 'tau' that we haven't stored coefficients in yet
  305. */
  306. c = gsl_matrix_subcolumn(A, i, i + 1, N - i - 1);
  307. hv = gsl_vector_subvector(tau, i + 1, N - (i + 1));
  308. gsl_vector_memcpy(&hv.vector, &c.vector);
  309. /* compute householder transformation of A(i+1:n,i) */
  310. tau_i = gsl_linalg_householder_transform(&hv.vector);
  311. /*
  312. * apply left householder matrix (I - tau_i v v') to
  313. * [ A | M_{23} ]
  314. */
  315. m = gsl_matrix_submatrix(M,
  316. top + i + 1,
  317. top + i,
  318. N - (i + 1),
  319. N_M - top - i);
  320. gsl_linalg_householder_hm(tau_i, &hv.vector, &m.matrix);
  321. /*
  322. * apply right householder matrix (I - tau_i v v') to
  323. *
  324. * [ M_{12} ]
  325. * [ A ]
  326. */
  327. m = gsl_matrix_submatrix(M,
  328. 0,
  329. top + i + 1,
  330. top + N,
  331. N - (i + 1));
  332. gsl_linalg_householder_mh(tau_i, &hv.vector, &m.matrix);
  333. /* save Householder coefficient */
  334. gsl_vector_set(tau, i, tau_i);
  335. /*
  336. * store Householder vector below the subdiagonal in column
  337. * i of the matrix. hv(1) does not need to be stored since
  338. * it is always 1.
  339. */
  340. c = gsl_vector_subvector(&c.vector, 1, c.vector.size - 1);
  341. hv = gsl_vector_subvector(&hv.vector, 1, hv.vector.size - 1);
  342. gsl_vector_memcpy(&c.vector, &hv.vector);
  343. }
  344. return GSL_SUCCESS;
  345. }
  346. } /* gsl_linalg_hessenberg_submatrix() */