gsl_eigen__jacobi.c 5.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260
  1. /* eigen/jacobi.c
  2. *
  3. * Copyright (C) 2004, 2007 Brian Gough, Gerard Jungman
  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_eigen.h"
  25. /* Algorithm 8.4.3 - Cyclic Jacobi. Golub & Van Loan, Matrix Computations */
  26. static inline double
  27. symschur2 (gsl_matrix * A, size_t p, size_t q, double *c, double *s)
  28. {
  29. double Apq = gsl_matrix_get (A, p, q);
  30. if (Apq != 0.0)
  31. {
  32. double App = gsl_matrix_get (A, p, p);
  33. double Aqq = gsl_matrix_get (A, q, q);
  34. double tau = (Aqq - App) / (2.0 * Apq);
  35. double t, c1;
  36. if (tau >= 0.0)
  37. {
  38. t = 1.0 / (tau + hypot (1.0, tau));
  39. }
  40. else
  41. {
  42. t = -1.0 / (-tau + hypot (1.0, tau));
  43. }
  44. c1 = 1.0 / hypot (1.0, t);
  45. *c = c1;
  46. *s = t * c1;
  47. }
  48. else
  49. {
  50. *c = 1.0;
  51. *s = 0.0;
  52. }
  53. /* reduction in off(A) is 2*(A_pq)^2 */
  54. return fabs (Apq);
  55. }
  56. inline static void
  57. apply_jacobi_L (gsl_matrix * A, size_t p, size_t q, double c, double s)
  58. {
  59. size_t j;
  60. const size_t N = A->size2;
  61. /* Apply rotation to matrix A, A' = J^T A */
  62. for (j = 0; j < N; j++)
  63. {
  64. double Apj = gsl_matrix_get (A, p, j);
  65. double Aqj = gsl_matrix_get (A, q, j);
  66. gsl_matrix_set (A, p, j, Apj * c - Aqj * s);
  67. gsl_matrix_set (A, q, j, Apj * s + Aqj * c);
  68. }
  69. }
  70. inline static void
  71. apply_jacobi_R (gsl_matrix * A, size_t p, size_t q, double c, double s)
  72. {
  73. size_t i;
  74. const size_t M = A->size1;
  75. /* Apply rotation to matrix A, A' = A J */
  76. for (i = 0; i < M; i++)
  77. {
  78. double Aip = gsl_matrix_get (A, i, p);
  79. double Aiq = gsl_matrix_get (A, i, q);
  80. gsl_matrix_set (A, i, p, Aip * c - Aiq * s);
  81. gsl_matrix_set (A, i, q, Aip * s + Aiq * c);
  82. }
  83. }
  84. inline static double
  85. norm (gsl_matrix * A)
  86. {
  87. size_t i, j, M = A->size1, N = A->size2;
  88. double sum = 0.0, scale = 0.0, ssq = 1.0;
  89. for (i = 0; i < M; i++)
  90. {
  91. for (j = 0; j < N; j++)
  92. {
  93. double Aij = gsl_matrix_get (A, i, j);
  94. if (Aij != 0.0)
  95. {
  96. double ax = fabs (Aij);
  97. if (scale < ax)
  98. {
  99. ssq = 1.0 + ssq * (scale / ax) * (scale / ax);
  100. scale = ax;
  101. }
  102. else
  103. {
  104. ssq += (ax / scale) * (ax / scale);
  105. }
  106. }
  107. }
  108. }
  109. sum = scale * sqrt (ssq);
  110. return sum;
  111. }
  112. int
  113. gsl_eigen_jacobi (gsl_matrix * a,
  114. gsl_vector * eval,
  115. gsl_matrix * evec, unsigned int max_rot, unsigned int *nrot)
  116. {
  117. size_t i, p, q;
  118. const size_t M = a->size1, N = a->size2;
  119. double red, redsum = 0.0;
  120. if (M != N)
  121. {
  122. GSL_ERROR ("eigenproblem requires square matrix", GSL_ENOTSQR);
  123. }
  124. else if (M != evec->size1 || M != evec->size2)
  125. {
  126. GSL_ERROR ("eigenvector matrix must match input matrix", GSL_EBADLEN);
  127. }
  128. else if (M != eval->size)
  129. {
  130. GSL_ERROR ("eigenvalue vector must match input matrix", GSL_EBADLEN);
  131. }
  132. gsl_vector_set_zero (eval);
  133. gsl_matrix_set_identity (evec);
  134. for (i = 0; i < max_rot; i++)
  135. {
  136. double nrm = norm (a);
  137. if (nrm == 0.0)
  138. break;
  139. for (p = 0; p < N; p++)
  140. {
  141. for (q = p + 1; q < N; q++)
  142. {
  143. double c, s;
  144. red = symschur2 (a, p, q, &c, &s);
  145. redsum += red;
  146. /* Compute A <- J^T A J */
  147. apply_jacobi_L (a, p, q, c, s);
  148. apply_jacobi_R (a, p, q, c, s);
  149. /* Compute V <- V J */
  150. apply_jacobi_R (evec, p, q, c, s);
  151. }
  152. }
  153. }
  154. *nrot = i;
  155. for (p = 0; p < N; p++)
  156. {
  157. double ep = gsl_matrix_get (a, p, p);
  158. gsl_vector_set (eval, p, ep);
  159. }
  160. if (i == max_rot)
  161. {
  162. return GSL_EMAXITER;
  163. }
  164. return GSL_SUCCESS;
  165. }
  166. int
  167. gsl_eigen_invert_jacobi (const gsl_matrix * a,
  168. gsl_matrix * ainv, unsigned int max_rot)
  169. {
  170. if (a->size1 != a->size2 || ainv->size1 != ainv->size2)
  171. {
  172. GSL_ERROR("jacobi method requires square matrix", GSL_ENOTSQR);
  173. }
  174. else if (a->size1 != ainv->size2)
  175. {
  176. GSL_ERROR ("inverse matrix must match input matrix", GSL_EBADLEN);
  177. }
  178. {
  179. const size_t n = a->size2;
  180. size_t i,j,k;
  181. unsigned int nrot = 0;
  182. int status;
  183. gsl_vector * eval = gsl_vector_alloc(n);
  184. gsl_matrix * evec = gsl_matrix_alloc(n, n);
  185. gsl_matrix * tmp = gsl_matrix_alloc(n, n);
  186. gsl_matrix_memcpy (tmp, a);
  187. status = gsl_eigen_jacobi(tmp, eval, evec, max_rot, &nrot);
  188. for(i=0; i<n; i++)
  189. {
  190. for(j=0; j<n; j++)
  191. {
  192. double ainv_ij = 0.0;
  193. for(k = 0; k<n; k++)
  194. {
  195. double f = 1.0 / gsl_vector_get(eval, k);
  196. double vik = gsl_matrix_get (evec, i, k);
  197. double vjk = gsl_matrix_get (evec, j, k);
  198. ainv_ij += vik * vjk * f;
  199. }
  200. gsl_matrix_set (ainv, i, j, ainv_ij);
  201. }
  202. }
  203. gsl_vector_free(eval);
  204. gsl_matrix_free(evec);
  205. gsl_matrix_free(tmp);
  206. if (status)
  207. {
  208. return status;
  209. }
  210. else
  211. {
  212. return GSL_SUCCESS;
  213. }
  214. }
  215. }