gsl_eigen__gensymmv.c 5.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201
  1. /* eigen/gensymmv.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. /*
  28. * This module computes the eigenvalues and eigenvectors of a real
  29. * generalized symmetric-definite eigensystem A x = \lambda B x, where
  30. * A and B are symmetric, and B is positive-definite.
  31. */
  32. static void gensymmv_normalize_eigenvectors(gsl_matrix *evec);
  33. /*
  34. gsl_eigen_gensymmv_alloc()
  35. Allocate a workspace for solving the generalized symmetric-definite
  36. eigenvalue problem. The size of this workspace is O(4n).
  37. Inputs: n - size of matrices
  38. Return: pointer to workspace
  39. */
  40. gsl_eigen_gensymmv_workspace *
  41. gsl_eigen_gensymmv_alloc(const size_t n)
  42. {
  43. gsl_eigen_gensymmv_workspace *w;
  44. if (n == 0)
  45. {
  46. GSL_ERROR_NULL ("matrix dimension must be positive integer",
  47. GSL_EINVAL);
  48. }
  49. w = (gsl_eigen_gensymmv_workspace *)calloc (1, sizeof (gsl_eigen_gensymmv_workspace));
  50. if (w == 0)
  51. {
  52. GSL_ERROR_NULL ("failed to allocate space for workspace", GSL_ENOMEM);
  53. }
  54. w->size = n;
  55. w->symmv_workspace_p = gsl_eigen_symmv_alloc(n);
  56. if (!w->symmv_workspace_p)
  57. {
  58. gsl_eigen_gensymmv_free(w);
  59. GSL_ERROR_NULL("failed to allocate space for symmv workspace", GSL_ENOMEM);
  60. }
  61. return (w);
  62. } /* gsl_eigen_gensymmv_alloc() */
  63. /*
  64. gsl_eigen_gensymmv_free()
  65. Free workspace w
  66. */
  67. void
  68. gsl_eigen_gensymmv_free (gsl_eigen_gensymmv_workspace * w)
  69. {
  70. if (w->symmv_workspace_p)
  71. gsl_eigen_symmv_free(w->symmv_workspace_p);
  72. free(w);
  73. } /* gsl_eigen_gensymmv_free() */
  74. /*
  75. gsl_eigen_gensymmv()
  76. Solve the generalized symmetric-definite eigenvalue problem
  77. A x = \lambda B x
  78. for the eigenvalues \lambda and eigenvectors x.
  79. Inputs: A - real symmetric matrix
  80. B - real symmetric and positive definite matrix
  81. eval - where to store eigenvalues
  82. evec - where to store eigenvectors
  83. w - workspace
  84. Return: success or error
  85. */
  86. int
  87. gsl_eigen_gensymmv (gsl_matrix * A, gsl_matrix * B, gsl_vector * eval,
  88. gsl_matrix * evec, gsl_eigen_gensymmv_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 (evec->size1 != evec->size2)
  105. {
  106. GSL_ERROR ("eigenvector matrix must be square", GSL_ENOTSQR);
  107. }
  108. else if (evec->size1 != N)
  109. {
  110. GSL_ERROR ("eigenvector matrix has wrong size", GSL_EBADLEN);
  111. }
  112. else if (w->size != N)
  113. {
  114. GSL_ERROR ("matrix size does not match workspace", GSL_EBADLEN);
  115. }
  116. else
  117. {
  118. int s;
  119. /* compute Cholesky factorization of B */
  120. s = gsl_linalg_cholesky_decomp(B);
  121. if (s != GSL_SUCCESS)
  122. return s; /* B is not positive definite */
  123. /* transform to standard symmetric eigenvalue problem */
  124. gsl_eigen_gensymm_standardize(A, B);
  125. /* compute eigenvalues and eigenvectors */
  126. s = gsl_eigen_symmv(A, eval, evec, w->symmv_workspace_p);
  127. if (s != GSL_SUCCESS)
  128. return s;
  129. /* backtransform eigenvectors: evec -> L^{-T} evec */
  130. gsl_blas_dtrsm(CblasLeft,
  131. CblasLower,
  132. CblasTrans,
  133. CblasNonUnit,
  134. 1.0,
  135. B,
  136. evec);
  137. /* the blas call destroyed the normalization - renormalize */
  138. gensymmv_normalize_eigenvectors(evec);
  139. return GSL_SUCCESS;
  140. }
  141. } /* gsl_eigen_gensymmv() */
  142. /********************************************
  143. * INTERNAL ROUTINES *
  144. ********************************************/
  145. /*
  146. gensymmv_normalize_eigenvectors()
  147. Normalize eigenvectors so that their Euclidean norm is 1
  148. Inputs: evec - eigenvectors
  149. */
  150. static void
  151. gensymmv_normalize_eigenvectors(gsl_matrix *evec)
  152. {
  153. const size_t N = evec->size1;
  154. size_t i; /* looping */
  155. for (i = 0; i < N; ++i)
  156. {
  157. gsl_vector_view vi = gsl_matrix_column(evec, i);
  158. double scale = 1.0 / gsl_blas_dnrm2(&vi.vector);
  159. gsl_blas_dscal(scale, &vi.vector);
  160. }
  161. } /* gensymmv_normalize_eigenvectors() */