gsl_eigen__hermv.c 6.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250
  1. /* eigen/hermv.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. #include "gsl__config.h"
  20. #include <stdlib.h>
  21. #include "gsl_math.h"
  22. #include "gsl_vector.h"
  23. #include "gsl_matrix.h"
  24. #include "gsl_complex_math.h"
  25. #include "gsl_linalg.h"
  26. #include "gsl_eigen.h"
  27. /* Compute eigenvalues/eigenvectors of complex hermitian matrix using
  28. reduction to real symmetric tridiagonal form, followed by QR
  29. iteration with implicit shifts.
  30. See Golub & Van Loan, "Matrix Computations" (3rd ed), Section 8.3 */
  31. #include "gsl_eigen__qrstep.c"
  32. gsl_eigen_hermv_workspace *
  33. gsl_eigen_hermv_alloc (const size_t n)
  34. {
  35. gsl_eigen_hermv_workspace * w ;
  36. if (n == 0)
  37. {
  38. GSL_ERROR_NULL ("matrix dimension must be positive integer", GSL_EINVAL);
  39. }
  40. w = (gsl_eigen_hermv_workspace *) malloc (sizeof(gsl_eigen_hermv_workspace));
  41. if (w == 0)
  42. {
  43. GSL_ERROR_NULL ("failed to allocate space for workspace", GSL_ENOMEM);
  44. }
  45. w->d = (double *) malloc (n * sizeof (double));
  46. if (w->d == 0)
  47. {
  48. free (w);
  49. GSL_ERROR_NULL ("failed to allocate space for diagonal", GSL_ENOMEM);
  50. }
  51. w->sd = (double *) malloc (n * sizeof (double));
  52. if (w->sd == 0)
  53. {
  54. free (w->d);
  55. free (w);
  56. GSL_ERROR_NULL ("failed to allocate space for subdiagonal", GSL_ENOMEM);
  57. }
  58. w->tau = (double *) malloc (2 * n * sizeof (double));
  59. if (w->tau == 0)
  60. {
  61. free (w->sd);
  62. free (w->d);
  63. free (w);
  64. GSL_ERROR_NULL ("failed to allocate space for tau", GSL_ENOMEM);
  65. }
  66. w->gc = (double *) malloc (n * sizeof (double));
  67. if (w->gc == 0)
  68. {
  69. free (w->tau);
  70. free (w->sd);
  71. free (w->d);
  72. free (w);
  73. GSL_ERROR_NULL ("failed to allocate space for cosines", GSL_ENOMEM);
  74. }
  75. w->gs = (double *) malloc (n * sizeof (double));
  76. if (w->gs == 0)
  77. {
  78. free (w->gc);
  79. free (w->tau);
  80. free (w->sd);
  81. free (w->d);
  82. free (w);
  83. GSL_ERROR_NULL ("failed to allocate space for sines", GSL_ENOMEM);
  84. }
  85. w->size = n;
  86. return w;
  87. }
  88. void
  89. gsl_eigen_hermv_free (gsl_eigen_hermv_workspace * w)
  90. {
  91. free (w->gs);
  92. free (w->gc);
  93. free (w->tau);
  94. free (w->sd);
  95. free (w->d);
  96. free (w);
  97. }
  98. int
  99. gsl_eigen_hermv (gsl_matrix_complex * A, gsl_vector * eval,
  100. gsl_matrix_complex * evec,
  101. gsl_eigen_hermv_workspace * w)
  102. {
  103. if (A->size1 != A->size2)
  104. {
  105. GSL_ERROR ("matrix must be square to compute eigenvalues", GSL_ENOTSQR);
  106. }
  107. else if (eval->size != A->size1)
  108. {
  109. GSL_ERROR ("eigenvalue vector must match matrix size", GSL_EBADLEN);
  110. }
  111. else if (evec->size1 != A->size1 || evec->size2 != A->size1)
  112. {
  113. GSL_ERROR ("eigenvector matrix must match matrix size", GSL_EBADLEN);
  114. }
  115. else
  116. {
  117. const size_t N = A->size1;
  118. double *const d = w->d;
  119. double *const sd = w->sd;
  120. size_t a, b;
  121. /* handle special case */
  122. if (N == 1)
  123. {
  124. gsl_complex A00 = gsl_matrix_complex_get (A, 0, 0);
  125. gsl_vector_set (eval, 0, GSL_REAL(A00));
  126. gsl_matrix_complex_set (evec, 0, 0, GSL_COMPLEX_ONE);
  127. return GSL_SUCCESS;
  128. }
  129. /* Transform the matrix into a symmetric tridiagonal form */
  130. {
  131. gsl_vector_view d_vec = gsl_vector_view_array (d, N);
  132. gsl_vector_view sd_vec = gsl_vector_view_array (sd, N - 1);
  133. gsl_vector_complex_view tau_vec = gsl_vector_complex_view_array (w->tau, N-1);
  134. gsl_linalg_hermtd_decomp (A, &tau_vec.vector);
  135. gsl_linalg_hermtd_unpack (A, &tau_vec.vector, evec, &d_vec.vector, &sd_vec.vector);
  136. }
  137. /* Make an initial pass through the tridiagonal decomposition
  138. to remove off-diagonal elements which are effectively zero */
  139. chop_small_elements (N, d, sd);
  140. /* Progressively reduce the matrix until it is diagonal */
  141. b = N - 1;
  142. while (b > 0)
  143. {
  144. if (sd[b - 1] == 0.0 || isnan(sd[b - 1]))
  145. {
  146. b--;
  147. continue;
  148. }
  149. /* Find the largest unreduced block (a,b) starting from b
  150. and working backwards */
  151. a = b - 1;
  152. while (a > 0)
  153. {
  154. if (sd[a - 1] == 0.0)
  155. {
  156. break;
  157. }
  158. a--;
  159. }
  160. {
  161. size_t i;
  162. const size_t n_block = b - a + 1;
  163. double *d_block = d + a;
  164. double *sd_block = sd + a;
  165. double * const gc = w->gc;
  166. double * const gs = w->gs;
  167. /* apply QR reduction with implicit deflation to the
  168. unreduced block */
  169. qrstep (n_block, d_block, sd_block, gc, gs);
  170. /* Apply Givens rotation Gij(c,s) to matrix Q, Q <- Q G */
  171. for (i = 0; i < n_block - 1; i++)
  172. {
  173. const double c = gc[i], s = gs[i];
  174. size_t k;
  175. for (k = 0; k < N; k++)
  176. {
  177. gsl_complex qki = gsl_matrix_complex_get (evec, k, a + i);
  178. gsl_complex qkj = gsl_matrix_complex_get (evec, k, a + i + 1);
  179. /* qki <= qki * c - qkj * s */
  180. /* qkj <= qki * s + qkj * c */
  181. gsl_complex x1 = gsl_complex_mul_real(qki, c);
  182. gsl_complex y1 = gsl_complex_mul_real(qkj, -s);
  183. gsl_complex x2 = gsl_complex_mul_real(qki, s);
  184. gsl_complex y2 = gsl_complex_mul_real(qkj, c);
  185. gsl_complex qqki = gsl_complex_add(x1, y1);
  186. gsl_complex qqkj = gsl_complex_add(x2, y2);
  187. gsl_matrix_complex_set (evec, k, a + i, qqki);
  188. gsl_matrix_complex_set (evec, k, a + i + 1, qqkj);
  189. }
  190. }
  191. /* remove any small off-diagonal elements */
  192. chop_small_elements (n_block, d_block, sd_block);
  193. }
  194. }
  195. {
  196. gsl_vector_view d_vec = gsl_vector_view_array (d, N);
  197. gsl_vector_memcpy (eval, &d_vec.vector);
  198. }
  199. return GSL_SUCCESS;
  200. }
  201. }