gsl_eigen__genherm.c 5.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226
  1. /* eigen/genherm.c
  2. *
  3. * Copyright (C) 2007 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. #include <stdlib.h>
  20. #include "gsl__config.h"
  21. #include "gsl_eigen.h"
  22. #include "gsl_linalg.h"
  23. #include "gsl_math.h"
  24. #include "gsl_blas.h"
  25. #include "gsl_vector.h"
  26. #include "gsl_matrix.h"
  27. #include "gsl_complex.h"
  28. #include "gsl_complex_math.h"
  29. /*
  30. * This module computes the eigenvalues of a complex generalized
  31. * hermitian-definite eigensystem A x = \lambda B x, where A and
  32. * B are hermitian, and B is positive-definite.
  33. */
  34. /*
  35. gsl_eigen_genherm_alloc()
  36. Allocate a workspace for solving the generalized hermitian-definite
  37. eigenvalue problem. The size of this workspace is O(3n).
  38. Inputs: n - size of matrices
  39. Return: pointer to workspace
  40. */
  41. gsl_eigen_genherm_workspace *
  42. gsl_eigen_genherm_alloc(const size_t n)
  43. {
  44. gsl_eigen_genherm_workspace *w;
  45. if (n == 0)
  46. {
  47. GSL_ERROR_NULL ("matrix dimension must be positive integer",
  48. GSL_EINVAL);
  49. }
  50. w = (gsl_eigen_genherm_workspace *) calloc (1, sizeof (gsl_eigen_genherm_workspace));
  51. if (w == 0)
  52. {
  53. GSL_ERROR_NULL ("failed to allocate space for workspace", GSL_ENOMEM);
  54. }
  55. w->size = n;
  56. w->herm_workspace_p = gsl_eigen_herm_alloc(n);
  57. if (!w->herm_workspace_p)
  58. {
  59. gsl_eigen_genherm_free(w);
  60. GSL_ERROR_NULL("failed to allocate space for herm workspace", GSL_ENOMEM);
  61. }
  62. return (w);
  63. } /* gsl_eigen_genherm_alloc() */
  64. /*
  65. gsl_eigen_genherm_free()
  66. Free workspace w
  67. */
  68. void
  69. gsl_eigen_genherm_free (gsl_eigen_genherm_workspace * w)
  70. {
  71. if (w->herm_workspace_p)
  72. gsl_eigen_herm_free(w->herm_workspace_p);
  73. free(w);
  74. } /* gsl_eigen_genherm_free() */
  75. /*
  76. gsl_eigen_genherm()
  77. Solve the generalized hermitian-definite eigenvalue problem
  78. A x = \lambda B x
  79. for the eigenvalues \lambda.
  80. Inputs: A - complex hermitian matrix
  81. B - complex hermitian and positive definite matrix
  82. eval - where to store eigenvalues
  83. w - workspace
  84. Return: success or error
  85. */
  86. int
  87. gsl_eigen_genherm (gsl_matrix_complex * A, gsl_matrix_complex * B,
  88. gsl_vector * eval, gsl_eigen_genherm_workspace * w)
  89. {
  90. const size_t N = A->size1;
  91. /* check matrix and vector sizes */
  92. if (N != A->size2)
  93. {
  94. GSL_ERROR ("matrix must be square to compute eigenvalues", GSL_ENOTSQR);
  95. }
  96. else if ((N != B->size1) || (N != B->size2))
  97. {
  98. GSL_ERROR ("B matrix dimensions must match A", GSL_EBADLEN);
  99. }
  100. else if (eval->size != N)
  101. {
  102. GSL_ERROR ("eigenvalue vector must match matrix size", GSL_EBADLEN);
  103. }
  104. else if (w->size != N)
  105. {
  106. GSL_ERROR ("matrix size does not match workspace", GSL_EBADLEN);
  107. }
  108. else
  109. {
  110. int s;
  111. /* compute Cholesky factorization of B */
  112. s = gsl_linalg_complex_cholesky_decomp(B);
  113. if (s != GSL_SUCCESS)
  114. return s; /* B is not positive definite */
  115. /* transform to standard hermitian eigenvalue problem */
  116. gsl_eigen_genherm_standardize(A, B);
  117. s = gsl_eigen_herm(A, eval, w->herm_workspace_p);
  118. return s;
  119. }
  120. } /* gsl_eigen_genherm() */
  121. /*
  122. gsl_eigen_genherm_standardize()
  123. Reduce the generalized hermitian-definite eigenproblem to
  124. the standard hermitian eigenproblem by computing
  125. C = L^{-1} A L^{-H}
  126. where L L^H is the Cholesky decomposition of B
  127. Inputs: A - (input/output) complex hermitian matrix
  128. B - complex hermitian, positive definite matrix in Cholesky form
  129. Return: success
  130. Notes: A is overwritten by L^{-1} A L^{-H}
  131. */
  132. int
  133. gsl_eigen_genherm_standardize(gsl_matrix_complex *A,
  134. const gsl_matrix_complex *B)
  135. {
  136. const size_t N = A->size1;
  137. size_t i;
  138. double a, b;
  139. gsl_complex y, z;
  140. GSL_SET_IMAG(&z, 0.0);
  141. for (i = 0; i < N; ++i)
  142. {
  143. /* update lower triangle of A(i:n, i:n) */
  144. y = gsl_matrix_complex_get(A, i, i);
  145. a = GSL_REAL(y);
  146. y = gsl_matrix_complex_get(B, i, i);
  147. b = GSL_REAL(y);
  148. a /= b * b;
  149. GSL_SET_REAL(&z, a);
  150. gsl_matrix_complex_set(A, i, i, z);
  151. if (i < N - 1)
  152. {
  153. gsl_vector_complex_view ai =
  154. gsl_matrix_complex_subcolumn(A, i, i + 1, N - i - 1);
  155. gsl_matrix_complex_view ma =
  156. gsl_matrix_complex_submatrix(A, i + 1, i + 1, N - i - 1, N - i - 1);
  157. gsl_vector_complex_const_view bi =
  158. gsl_matrix_complex_const_subcolumn(B, i, i + 1, N - i - 1);
  159. gsl_matrix_complex_const_view mb =
  160. gsl_matrix_complex_const_submatrix(B, i + 1, i + 1, N - i - 1, N - i - 1);
  161. gsl_blas_zdscal(1.0 / b, &ai.vector);
  162. GSL_SET_REAL(&z, -0.5 * a);
  163. gsl_blas_zaxpy(z, &bi.vector, &ai.vector);
  164. gsl_blas_zher2(CblasLower,
  165. GSL_COMPLEX_NEGONE,
  166. &ai.vector,
  167. &bi.vector,
  168. &ma.matrix);
  169. gsl_blas_zaxpy(z, &bi.vector, &ai.vector);
  170. gsl_blas_ztrsv(CblasLower,
  171. CblasNoTrans,
  172. CblasNonUnit,
  173. &mb.matrix,
  174. &ai.vector);
  175. }
  176. }
  177. return GSL_SUCCESS;
  178. } /* gsl_eigen_genherm_standardize() */