gsl_eigen__symmv.c 5.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216
  1. /* eigen/symmv.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_linalg.h"
  25. #include "gsl_eigen.h"
  26. /* Compute eigenvalues/eigenvectors of real symmetric matrix using
  27. reduction to tridiagonal form, followed by QR iteration with
  28. implicit shifts.
  29. See Golub & Van Loan, "Matrix Computations" (3rd ed), Section 8.3
  30. */
  31. #include "gsl_eigen__qrstep.c"
  32. gsl_eigen_symmv_workspace *
  33. gsl_eigen_symmv_alloc (const size_t n)
  34. {
  35. gsl_eigen_symmv_workspace * w ;
  36. if (n == 0)
  37. {
  38. GSL_ERROR_NULL ("matrix dimension must be positive integer", GSL_EINVAL);
  39. }
  40. w= ((gsl_eigen_symmv_workspace *) malloc (sizeof(gsl_eigen_symmv_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. GSL_ERROR_NULL ("failed to allocate space for diagonal", GSL_ENOMEM);
  49. }
  50. w->sd = (double *) malloc (n * sizeof (double));
  51. if (w->sd == 0)
  52. {
  53. GSL_ERROR_NULL ("failed to allocate space for subdiagonal", GSL_ENOMEM);
  54. }
  55. w->gc = (double *) malloc (n * sizeof (double));
  56. if (w->gc == 0)
  57. {
  58. GSL_ERROR_NULL ("failed to allocate space for cosines", GSL_ENOMEM);
  59. }
  60. w->gs = (double *) malloc (n * sizeof (double));
  61. if (w->gs == 0)
  62. {
  63. GSL_ERROR_NULL ("failed to allocate space for sines", GSL_ENOMEM);
  64. }
  65. w->size = n;
  66. return w;
  67. }
  68. void
  69. gsl_eigen_symmv_free (gsl_eigen_symmv_workspace * w)
  70. {
  71. free(w->gs);
  72. free(w->gc);
  73. free(w->sd);
  74. free(w->d);
  75. free(w);
  76. }
  77. int
  78. gsl_eigen_symmv (gsl_matrix * A, gsl_vector * eval, gsl_matrix * evec,
  79. gsl_eigen_symmv_workspace * w)
  80. {
  81. if (A->size1 != A->size2)
  82. {
  83. GSL_ERROR ("matrix must be square to compute eigenvalues", GSL_ENOTSQR);
  84. }
  85. else if (eval->size != A->size1)
  86. {
  87. GSL_ERROR ("eigenvalue vector must match matrix size", GSL_EBADLEN);
  88. }
  89. else if (evec->size1 != A->size1 || evec->size2 != A->size1)
  90. {
  91. GSL_ERROR ("eigenvector matrix must match matrix size", GSL_EBADLEN);
  92. }
  93. else
  94. {
  95. double *const d = w->d;
  96. double *const sd = w->sd;
  97. const size_t N = A->size1;
  98. size_t a, b;
  99. /* handle special case */
  100. if (N == 1)
  101. {
  102. double A00 = gsl_matrix_get (A, 0, 0);
  103. gsl_vector_set (eval, 0, A00);
  104. gsl_matrix_set (evec, 0, 0, 1.0);
  105. return GSL_SUCCESS;
  106. }
  107. /* use sd as the temporary workspace for the decomposition when
  108. computing eigenvectors */
  109. {
  110. gsl_vector_view d_vec = gsl_vector_view_array (d, N);
  111. gsl_vector_view sd_vec = gsl_vector_view_array (sd, N - 1);
  112. gsl_vector_view tau = gsl_vector_view_array (sd, N - 1);
  113. gsl_linalg_symmtd_decomp (A, &tau.vector);
  114. gsl_linalg_symmtd_unpack (A, &tau.vector, evec, &d_vec.vector, &sd_vec.vector);
  115. }
  116. /* Make an initial pass through the tridiagonal decomposition
  117. to remove off-diagonal elements which are effectively zero */
  118. chop_small_elements (N, d, sd);
  119. /* Progressively reduce the matrix until it is diagonal */
  120. b = N - 1;
  121. while (b > 0)
  122. {
  123. if (sd[b - 1] == 0.0 || isnan(sd[b - 1]))
  124. {
  125. b--;
  126. continue;
  127. }
  128. /* Find the largest unreduced block (a,b) starting from b
  129. and working backwards */
  130. a = b - 1;
  131. while (a > 0)
  132. {
  133. if (sd[a - 1] == 0.0)
  134. {
  135. break;
  136. }
  137. a--;
  138. }
  139. {
  140. size_t i;
  141. const size_t n_block = b - a + 1;
  142. double *d_block = d + a;
  143. double *sd_block = sd + a;
  144. double * const gc = w->gc;
  145. double * const gs = w->gs;
  146. /* apply QR reduction with implicit deflation to the
  147. unreduced block */
  148. qrstep (n_block, d_block, sd_block, gc, gs);
  149. /* Apply Givens rotation Gij(c,s) to matrix Q, Q <- Q G */
  150. for (i = 0; i < n_block - 1; i++)
  151. {
  152. const double c = gc[i], s = gs[i];
  153. size_t k;
  154. for (k = 0; k < N; k++)
  155. {
  156. double qki = gsl_matrix_get (evec, k, a + i);
  157. double qkj = gsl_matrix_get (evec, k, a + i + 1);
  158. gsl_matrix_set (evec, k, a + i, qki * c - qkj * s);
  159. gsl_matrix_set (evec, k, a + i + 1, qki * s + qkj * c);
  160. }
  161. }
  162. /* remove any small off-diagonal elements */
  163. chop_small_elements (N, d, sd);
  164. }
  165. }
  166. {
  167. gsl_vector_view d_vec = gsl_vector_view_array (d, N);
  168. gsl_vector_memcpy (eval, &d_vec.vector);
  169. }
  170. return GSL_SUCCESS;
  171. }
  172. }