gsl_eigen__sort.c 9.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347
  1. /* eigen/sort.c
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2006, 2007 Gerard Jungman, 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. /* Author: G. Jungman, Modified: B. Gough. */
  20. #include "gsl__config.h"
  21. #include <stdlib.h>
  22. #include "gsl_math.h"
  23. #include "gsl_eigen.h"
  24. #include "gsl_complex.h"
  25. #include "gsl_complex_math.h"
  26. /* The eigen_sort below is not very good, but it is simple and
  27. * self-contained. We can always implement an improved sort later. */
  28. int
  29. gsl_eigen_symmv_sort (gsl_vector * eval, gsl_matrix * evec,
  30. gsl_eigen_sort_t sort_type)
  31. {
  32. if (evec->size1 != evec->size2)
  33. {
  34. GSL_ERROR ("eigenvector matrix must be square", GSL_ENOTSQR);
  35. }
  36. else if (eval->size != evec->size1)
  37. {
  38. GSL_ERROR ("eigenvalues must match eigenvector matrix", GSL_EBADLEN);
  39. }
  40. else
  41. {
  42. const size_t N = eval->size;
  43. size_t i;
  44. for (i = 0; i < N - 1; i++)
  45. {
  46. size_t j;
  47. size_t k = i;
  48. double ek = gsl_vector_get (eval, i);
  49. /* search for something to swap */
  50. for (j = i + 1; j < N; j++)
  51. {
  52. int test;
  53. const double ej = gsl_vector_get (eval, j);
  54. switch (sort_type)
  55. {
  56. case GSL_EIGEN_SORT_VAL_ASC:
  57. test = (ej < ek);
  58. break;
  59. case GSL_EIGEN_SORT_VAL_DESC:
  60. test = (ej > ek);
  61. break;
  62. case GSL_EIGEN_SORT_ABS_ASC:
  63. test = (fabs (ej) < fabs (ek));
  64. break;
  65. case GSL_EIGEN_SORT_ABS_DESC:
  66. test = (fabs (ej) > fabs (ek));
  67. break;
  68. default:
  69. GSL_ERROR ("unrecognized sort type", GSL_EINVAL);
  70. }
  71. if (test)
  72. {
  73. k = j;
  74. ek = ej;
  75. }
  76. }
  77. if (k != i)
  78. {
  79. /* swap eigenvalues */
  80. gsl_vector_swap_elements (eval, i, k);
  81. /* swap eigenvectors */
  82. gsl_matrix_swap_columns (evec, i, k);
  83. }
  84. }
  85. return GSL_SUCCESS;
  86. }
  87. }
  88. int
  89. gsl_eigen_hermv_sort (gsl_vector * eval, gsl_matrix_complex * evec,
  90. gsl_eigen_sort_t sort_type)
  91. {
  92. if (evec->size1 != evec->size2)
  93. {
  94. GSL_ERROR ("eigenvector matrix must be square", GSL_ENOTSQR);
  95. }
  96. else if (eval->size != evec->size1)
  97. {
  98. GSL_ERROR ("eigenvalues must match eigenvector matrix", GSL_EBADLEN);
  99. }
  100. else
  101. {
  102. const size_t N = eval->size;
  103. size_t i;
  104. for (i = 0; i < N - 1; i++)
  105. {
  106. size_t j;
  107. size_t k = i;
  108. double ek = gsl_vector_get (eval, i);
  109. /* search for something to swap */
  110. for (j = i + 1; j < N; j++)
  111. {
  112. int test;
  113. const double ej = gsl_vector_get (eval, j);
  114. switch (sort_type)
  115. {
  116. case GSL_EIGEN_SORT_VAL_ASC:
  117. test = (ej < ek);
  118. break;
  119. case GSL_EIGEN_SORT_VAL_DESC:
  120. test = (ej > ek);
  121. break;
  122. case GSL_EIGEN_SORT_ABS_ASC:
  123. test = (fabs (ej) < fabs (ek));
  124. break;
  125. case GSL_EIGEN_SORT_ABS_DESC:
  126. test = (fabs (ej) > fabs (ek));
  127. break;
  128. default:
  129. GSL_ERROR ("unrecognized sort type", GSL_EINVAL);
  130. }
  131. if (test)
  132. {
  133. k = j;
  134. ek = ej;
  135. }
  136. }
  137. if (k != i)
  138. {
  139. /* swap eigenvalues */
  140. gsl_vector_swap_elements (eval, i, k);
  141. /* swap eigenvectors */
  142. gsl_matrix_complex_swap_columns (evec, i, k);
  143. }
  144. }
  145. return GSL_SUCCESS;
  146. }
  147. }
  148. int
  149. gsl_eigen_nonsymmv_sort (gsl_vector_complex * eval,
  150. gsl_matrix_complex * evec,
  151. gsl_eigen_sort_t sort_type)
  152. {
  153. if (evec && (evec->size1 != evec->size2))
  154. {
  155. GSL_ERROR ("eigenvector matrix must be square", GSL_ENOTSQR);
  156. }
  157. else if (evec && (eval->size != evec->size1))
  158. {
  159. GSL_ERROR ("eigenvalues must match eigenvector matrix", GSL_EBADLEN);
  160. }
  161. else
  162. {
  163. const size_t N = eval->size;
  164. size_t i;
  165. for (i = 0; i < N - 1; i++)
  166. {
  167. size_t j;
  168. size_t k = i;
  169. gsl_complex ek = gsl_vector_complex_get (eval, i);
  170. /* search for something to swap */
  171. for (j = i + 1; j < N; j++)
  172. {
  173. int test;
  174. const gsl_complex ej = gsl_vector_complex_get (eval, j);
  175. switch (sort_type)
  176. {
  177. case GSL_EIGEN_SORT_ABS_ASC:
  178. test = (gsl_complex_abs (ej) < gsl_complex_abs (ek));
  179. break;
  180. case GSL_EIGEN_SORT_ABS_DESC:
  181. test = (gsl_complex_abs (ej) > gsl_complex_abs (ek));
  182. break;
  183. case GSL_EIGEN_SORT_VAL_ASC:
  184. case GSL_EIGEN_SORT_VAL_DESC:
  185. default:
  186. GSL_ERROR ("invalid sort type", GSL_EINVAL);
  187. }
  188. if (test)
  189. {
  190. k = j;
  191. ek = ej;
  192. }
  193. }
  194. if (k != i)
  195. {
  196. /* swap eigenvalues */
  197. gsl_vector_complex_swap_elements (eval, i, k);
  198. /* swap eigenvectors */
  199. if (evec)
  200. gsl_matrix_complex_swap_columns (evec, i, k);
  201. }
  202. }
  203. return GSL_SUCCESS;
  204. }
  205. }
  206. int
  207. gsl_eigen_gensymmv_sort (gsl_vector * eval, gsl_matrix * evec,
  208. gsl_eigen_sort_t sort_type)
  209. {
  210. int s;
  211. s = gsl_eigen_symmv_sort(eval, evec, sort_type);
  212. return s;
  213. }
  214. int
  215. gsl_eigen_genhermv_sort (gsl_vector * eval, gsl_matrix_complex * evec,
  216. gsl_eigen_sort_t sort_type)
  217. {
  218. int s;
  219. s = gsl_eigen_hermv_sort(eval, evec, sort_type);
  220. return s;
  221. }
  222. int
  223. gsl_eigen_genv_sort (gsl_vector_complex * alpha, gsl_vector * beta,
  224. gsl_matrix_complex * evec, gsl_eigen_sort_t sort_type)
  225. {
  226. if (evec->size1 != evec->size2)
  227. {
  228. GSL_ERROR ("eigenvector matrix must be square", GSL_ENOTSQR);
  229. }
  230. else if (alpha->size != evec->size1 || beta->size != evec->size1)
  231. {
  232. GSL_ERROR ("eigenvalues must match eigenvector matrix", GSL_EBADLEN);
  233. }
  234. else
  235. {
  236. const size_t N = alpha->size;
  237. size_t i;
  238. for (i = 0; i < N - 1; i++)
  239. {
  240. size_t j;
  241. size_t k = i;
  242. gsl_complex ak = gsl_vector_complex_get (alpha, i);
  243. double bk = gsl_vector_get(beta, i);
  244. gsl_complex ek;
  245. if (bk < GSL_DBL_EPSILON)
  246. {
  247. GSL_SET_COMPLEX(&ek,
  248. GSL_SIGN(GSL_REAL(ak)) ? GSL_POSINF : GSL_NEGINF,
  249. GSL_SIGN(GSL_IMAG(ak)) ? GSL_POSINF : GSL_NEGINF);
  250. }
  251. else
  252. ek = gsl_complex_div_real(ak, bk);
  253. /* search for something to swap */
  254. for (j = i + 1; j < N; j++)
  255. {
  256. int test;
  257. const gsl_complex aj = gsl_vector_complex_get (alpha, j);
  258. double bj = gsl_vector_get(beta, j);
  259. gsl_complex ej;
  260. if (bj < GSL_DBL_EPSILON)
  261. {
  262. GSL_SET_COMPLEX(&ej,
  263. GSL_SIGN(GSL_REAL(aj)) ? GSL_POSINF : GSL_NEGINF,
  264. GSL_SIGN(GSL_IMAG(aj)) ? GSL_POSINF : GSL_NEGINF);
  265. }
  266. else
  267. ej = gsl_complex_div_real(aj, bj);
  268. switch (sort_type)
  269. {
  270. case GSL_EIGEN_SORT_ABS_ASC:
  271. test = (gsl_complex_abs (ej) < gsl_complex_abs (ek));
  272. break;
  273. case GSL_EIGEN_SORT_ABS_DESC:
  274. test = (gsl_complex_abs (ej) > gsl_complex_abs (ek));
  275. break;
  276. case GSL_EIGEN_SORT_VAL_ASC:
  277. case GSL_EIGEN_SORT_VAL_DESC:
  278. default:
  279. GSL_ERROR ("invalid sort type", GSL_EINVAL);
  280. }
  281. if (test)
  282. {
  283. k = j;
  284. ek = ej;
  285. }
  286. }
  287. if (k != i)
  288. {
  289. /* swap eigenvalues */
  290. gsl_vector_complex_swap_elements (alpha, i, k);
  291. gsl_vector_swap_elements (beta, i, k);
  292. /* swap eigenvectors */
  293. gsl_matrix_complex_swap_columns (evec, i, k);
  294. }
  295. }
  296. return GSL_SUCCESS;
  297. }
  298. }