gsl_eigen__genhermv.c 5.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204
  1. /* eigen/genhermv.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 and eigenvectors of a complex
  31. * generalized hermitian-definite eigensystem A x = \lambda B x, where
  32. * A and B are hermitian, and B is positive-definite.
  33. */
  34. static void genhermv_normalize_eigenvectors(gsl_matrix_complex *evec);
  35. /*
  36. gsl_eigen_genhermv_alloc()
  37. Allocate a workspace for solving the generalized hermitian-definite
  38. eigenvalue problem. The size of this workspace is O(5n).
  39. Inputs: n - size of matrices
  40. Return: pointer to workspace
  41. */
  42. gsl_eigen_genhermv_workspace *
  43. gsl_eigen_genhermv_alloc(const size_t n)
  44. {
  45. gsl_eigen_genhermv_workspace *w;
  46. if (n == 0)
  47. {
  48. GSL_ERROR_NULL ("matrix dimension must be positive integer",
  49. GSL_EINVAL);
  50. }
  51. w = (gsl_eigen_genhermv_workspace *) calloc (1, sizeof (gsl_eigen_genhermv_workspace));
  52. if (w == 0)
  53. {
  54. GSL_ERROR_NULL ("failed to allocate space for workspace", GSL_ENOMEM);
  55. }
  56. w->size = n;
  57. w->hermv_workspace_p = gsl_eigen_hermv_alloc(n);
  58. if (!w->hermv_workspace_p)
  59. {
  60. gsl_eigen_genhermv_free(w);
  61. GSL_ERROR_NULL("failed to allocate space for hermv workspace", GSL_ENOMEM);
  62. }
  63. return (w);
  64. } /* gsl_eigen_genhermv_alloc() */
  65. /*
  66. gsl_eigen_genhermv_free()
  67. Free workspace w
  68. */
  69. void
  70. gsl_eigen_genhermv_free (gsl_eigen_genhermv_workspace * w)
  71. {
  72. if (w->hermv_workspace_p)
  73. gsl_eigen_hermv_free(w->hermv_workspace_p);
  74. free(w);
  75. } /* gsl_eigen_genhermv_free() */
  76. /*
  77. gsl_eigen_genhermv()
  78. Solve the generalized hermitian-definite eigenvalue problem
  79. A x = \lambda B x
  80. for the eigenvalues \lambda and eigenvectors x.
  81. Inputs: A - complex hermitian matrix
  82. B - complex hermitian and positive definite matrix
  83. eval - where to store eigenvalues
  84. evec - where to store eigenvectors
  85. w - workspace
  86. Return: success or error
  87. */
  88. int
  89. gsl_eigen_genhermv (gsl_matrix_complex * A, gsl_matrix_complex * B,
  90. gsl_vector * eval, gsl_matrix_complex * evec,
  91. gsl_eigen_genhermv_workspace * w)
  92. {
  93. const size_t N = A->size1;
  94. /* check matrix and vector sizes */
  95. if (N != A->size2)
  96. {
  97. GSL_ERROR ("matrix must be square to compute eigenvalues", GSL_ENOTSQR);
  98. }
  99. else if ((N != B->size1) || (N != B->size2))
  100. {
  101. GSL_ERROR ("B matrix dimensions must match A", GSL_EBADLEN);
  102. }
  103. else if (eval->size != N)
  104. {
  105. GSL_ERROR ("eigenvalue vector must match matrix size", GSL_EBADLEN);
  106. }
  107. else if (evec->size1 != evec->size2)
  108. {
  109. GSL_ERROR ("eigenvector matrix must be square", GSL_ENOTSQR);
  110. }
  111. else if (evec->size1 != N)
  112. {
  113. GSL_ERROR ("eigenvector matrix has wrong size", GSL_EBADLEN);
  114. }
  115. else if (w->size != N)
  116. {
  117. GSL_ERROR ("matrix size does not match workspace", GSL_EBADLEN);
  118. }
  119. else
  120. {
  121. int s;
  122. /* compute Cholesky factorization of B */
  123. s = gsl_linalg_complex_cholesky_decomp(B);
  124. if (s != GSL_SUCCESS)
  125. return s; /* B is not positive definite */
  126. /* transform to standard hermitian eigenvalue problem */
  127. gsl_eigen_genherm_standardize(A, B);
  128. /* compute eigenvalues and eigenvectors */
  129. s = gsl_eigen_hermv(A, eval, evec, w->hermv_workspace_p);
  130. if (s != GSL_SUCCESS)
  131. return s;
  132. /* backtransform eigenvectors: evec -> L^{-H} evec */
  133. gsl_blas_ztrsm(CblasLeft,
  134. CblasLower,
  135. CblasConjTrans,
  136. CblasNonUnit,
  137. GSL_COMPLEX_ONE,
  138. B,
  139. evec);
  140. /* the blas call destroyed the normalization - renormalize */
  141. genhermv_normalize_eigenvectors(evec);
  142. return GSL_SUCCESS;
  143. }
  144. } /* gsl_eigen_genhermv() */
  145. /********************************************
  146. * INTERNAL ROUTINES *
  147. ********************************************/
  148. /*
  149. genhermv_normalize_eigenvectors()
  150. Normalize eigenvectors so that their Euclidean norm is 1
  151. Inputs: evec - eigenvectors
  152. */
  153. static void
  154. genhermv_normalize_eigenvectors(gsl_matrix_complex *evec)
  155. {
  156. const size_t N = evec->size1;
  157. size_t i; /* looping */
  158. for (i = 0; i < N; ++i)
  159. {
  160. gsl_vector_complex_view vi = gsl_matrix_complex_column(evec, i);
  161. double scale = 1.0 / gsl_blas_dznrm2(&vi.vector);
  162. gsl_blas_zdscal(scale, &vi.vector);
  163. }
  164. } /* genhermv_normalize_eigenvectors() */