gsl_linalg__hermtd.c 7.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241
  1. /* linalg/hermtd.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 hermitian matrix A into
  20. *
  21. * A = U T U'
  22. *
  23. * where U is unitary and T is real symmetric tridiagonal. Only the
  24. * diagonal and lower triangular part of A is referenced and modified.
  25. *
  26. * On exit, T is stored in the diagonal and first subdiagonal of
  27. * A. Since T is symmetric the upper diagonal is not stored.
  28. *
  29. * U is stored as a packed set of Householder transformations in the
  30. * lower triangular part of the input matrix below the first subdiagonal.
  31. *
  32. * The full matrix for Q can be obtained as the product
  33. *
  34. * Q = Q_N ... Q_2 Q_1
  35. *
  36. * where
  37. *
  38. * Q_i = (I - tau_i * v_i * v_i')
  39. *
  40. * and where v_i is a Householder vector
  41. *
  42. * v_i = [0, ..., 0, 1, A(i+2,i), A(i+3,i), ... , A(N,i)]
  43. *
  44. * This storage scheme is the same as in LAPACK. See LAPACK's
  45. * chetd2.f for details.
  46. *
  47. * See Golub & Van Loan, "Matrix Computations" (3rd ed), Section 8.3 */
  48. #include "gsl__config.h"
  49. #include <stdlib.h>
  50. #include "gsl_math.h"
  51. #include "gsl_vector.h"
  52. #include "gsl_matrix.h"
  53. #include "gsl_blas.h"
  54. #include "gsl_complex_math.h"
  55. #include "gsl_linalg.h"
  56. int
  57. gsl_linalg_hermtd_decomp (gsl_matrix_complex * A, gsl_vector_complex * tau)
  58. {
  59. if (A->size1 != A->size2)
  60. {
  61. GSL_ERROR ("hermitian tridiagonal decomposition requires square matrix",
  62. GSL_ENOTSQR);
  63. }
  64. else if (tau->size + 1 != A->size1)
  65. {
  66. GSL_ERROR ("size of tau must be (matrix size - 1)", GSL_EBADLEN);
  67. }
  68. else
  69. {
  70. const size_t N = A->size1;
  71. size_t i;
  72. const gsl_complex zero = gsl_complex_rect (0.0, 0.0);
  73. const gsl_complex one = gsl_complex_rect (1.0, 0.0);
  74. const gsl_complex neg_one = gsl_complex_rect (-1.0, 0.0);
  75. for (i = 0 ; i < N - 1; i++)
  76. {
  77. gsl_vector_complex_view c = gsl_matrix_complex_column (A, i);
  78. gsl_vector_complex_view v = gsl_vector_complex_subvector (&c.vector, i + 1, N - (i + 1));
  79. gsl_complex tau_i = gsl_linalg_complex_householder_transform (&v.vector);
  80. /* Apply the transformation H^T A H to the remaining columns */
  81. if ((i + 1) < (N - 1)
  82. && !(GSL_REAL(tau_i) == 0.0 && GSL_IMAG(tau_i) == 0.0))
  83. {
  84. gsl_matrix_complex_view m =
  85. gsl_matrix_complex_submatrix (A, i + 1, i + 1,
  86. N - (i+1), N - (i+1));
  87. gsl_complex ei = gsl_vector_complex_get(&v.vector, 0);
  88. gsl_vector_complex_view x = gsl_vector_complex_subvector (tau, i, N-(i+1));
  89. gsl_vector_complex_set (&v.vector, 0, one);
  90. /* x = tau * A * v */
  91. gsl_blas_zhemv (CblasLower, tau_i, &m.matrix, &v.vector, zero, &x.vector);
  92. /* w = x - (1/2) tau * (x' * v) * v */
  93. {
  94. gsl_complex xv, txv, alpha;
  95. gsl_blas_zdotc(&x.vector, &v.vector, &xv);
  96. txv = gsl_complex_mul(tau_i, xv);
  97. alpha = gsl_complex_mul_real(txv, -0.5);
  98. gsl_blas_zaxpy(alpha, &v.vector, &x.vector);
  99. }
  100. /* apply the transformation A = A - v w' - w v' */
  101. gsl_blas_zher2(CblasLower, neg_one, &v.vector, &x.vector, &m.matrix);
  102. gsl_vector_complex_set (&v.vector, 0, ei);
  103. }
  104. gsl_vector_complex_set (tau, i, tau_i);
  105. }
  106. return GSL_SUCCESS;
  107. }
  108. }
  109. /* Form the orthogonal matrix Q from the packed QR matrix */
  110. int
  111. gsl_linalg_hermtd_unpack (const gsl_matrix_complex * A,
  112. const gsl_vector_complex * tau,
  113. gsl_matrix_complex * Q,
  114. gsl_vector * diag,
  115. gsl_vector * sdiag)
  116. {
  117. if (A->size1 != A->size2)
  118. {
  119. GSL_ERROR ("matrix A must be sqaure", GSL_ENOTSQR);
  120. }
  121. else if (tau->size + 1 != A->size1)
  122. {
  123. GSL_ERROR ("size of tau must be (matrix size - 1)", GSL_EBADLEN);
  124. }
  125. else if (Q->size1 != A->size1 || Q->size2 != A->size1)
  126. {
  127. GSL_ERROR ("size of Q must match size of A", GSL_EBADLEN);
  128. }
  129. else if (diag->size != A->size1)
  130. {
  131. GSL_ERROR ("size of diagonal must match size of A", GSL_EBADLEN);
  132. }
  133. else if (sdiag->size + 1 != A->size1)
  134. {
  135. GSL_ERROR ("size of subdiagonal must be (matrix size - 1)", GSL_EBADLEN);
  136. }
  137. else
  138. {
  139. const size_t N = A->size1;
  140. size_t i;
  141. /* Initialize Q to the identity */
  142. gsl_matrix_complex_set_identity (Q);
  143. for (i = N - 1; i > 0 && i--;)
  144. {
  145. gsl_complex ti = gsl_vector_complex_get (tau, i);
  146. gsl_vector_complex_const_view c = gsl_matrix_complex_const_column (A, i);
  147. gsl_vector_complex_const_view h =
  148. gsl_vector_complex_const_subvector (&c.vector, i + 1, N - (i+1));
  149. gsl_matrix_complex_view m =
  150. gsl_matrix_complex_submatrix (Q, i + 1, i + 1, N-(i+1), N-(i+1));
  151. gsl_linalg_complex_householder_hm (ti, &h.vector, &m.matrix);
  152. }
  153. /* Copy diagonal into diag */
  154. for (i = 0; i < N; i++)
  155. {
  156. gsl_complex Aii = gsl_matrix_complex_get (A, i, i);
  157. gsl_vector_set (diag, i, GSL_REAL(Aii));
  158. }
  159. /* Copy subdiagonal into sdiag */
  160. for (i = 0; i < N - 1; i++)
  161. {
  162. gsl_complex Aji = gsl_matrix_complex_get (A, i+1, i);
  163. gsl_vector_set (sdiag, i, GSL_REAL(Aji));
  164. }
  165. return GSL_SUCCESS;
  166. }
  167. }
  168. int
  169. gsl_linalg_hermtd_unpack_T (const gsl_matrix_complex * A,
  170. gsl_vector * diag,
  171. gsl_vector * sdiag)
  172. {
  173. if (A->size1 != A->size2)
  174. {
  175. GSL_ERROR ("matrix A must be sqaure", GSL_ENOTSQR);
  176. }
  177. else if (diag->size != A->size1)
  178. {
  179. GSL_ERROR ("size of diagonal must match size of A", GSL_EBADLEN);
  180. }
  181. else if (sdiag->size + 1 != A->size1)
  182. {
  183. GSL_ERROR ("size of subdiagonal must be (matrix size - 1)", GSL_EBADLEN);
  184. }
  185. else
  186. {
  187. const size_t N = A->size1;
  188. size_t i;
  189. /* Copy diagonal into diag */
  190. for (i = 0; i < N; i++)
  191. {
  192. gsl_complex Aii = gsl_matrix_complex_get (A, i, i);
  193. gsl_vector_set (diag, i, GSL_REAL(Aii));
  194. }
  195. /* Copy subdiagonal into sd */
  196. for (i = 0; i < N - 1; i++)
  197. {
  198. gsl_complex Aji = gsl_matrix_complex_get (A, i+1, i);
  199. gsl_vector_set (sdiag, i, GSL_REAL(Aji));
  200. }
  201. return GSL_SUCCESS;
  202. }
  203. }