gsl_eigen__gensymm.c 5.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212
  1. /* eigen/gensymm.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 of a real generalized
  29. * symmetric-definite eigensystem A x = \lambda B x, where A and
  30. * B are symmetric, and B is positive-definite.
  31. */
  32. /*
  33. gsl_eigen_gensymm_alloc()
  34. Allocate a workspace for solving the generalized symmetric-definite
  35. eigenvalue problem. The size of this workspace is O(2n).
  36. Inputs: n - size of matrices
  37. Return: pointer to workspace
  38. */
  39. gsl_eigen_gensymm_workspace *
  40. gsl_eigen_gensymm_alloc(const size_t n)
  41. {
  42. gsl_eigen_gensymm_workspace *w;
  43. if (n == 0)
  44. {
  45. GSL_ERROR_NULL ("matrix dimension must be positive integer",
  46. GSL_EINVAL);
  47. }
  48. w = (gsl_eigen_gensymm_workspace *) calloc (1, sizeof (gsl_eigen_gensymm_workspace));
  49. if (w == 0)
  50. {
  51. GSL_ERROR_NULL ("failed to allocate space for workspace", GSL_ENOMEM);
  52. }
  53. w->size = n;
  54. w->symm_workspace_p = gsl_eigen_symm_alloc(n);
  55. if (!w->symm_workspace_p)
  56. {
  57. gsl_eigen_gensymm_free(w);
  58. GSL_ERROR_NULL("failed to allocate space for symm workspace", GSL_ENOMEM);
  59. }
  60. return (w);
  61. } /* gsl_eigen_gensymm_alloc() */
  62. /*
  63. gsl_eigen_gensymm_free()
  64. Free workspace w
  65. */
  66. void
  67. gsl_eigen_gensymm_free (gsl_eigen_gensymm_workspace * w)
  68. {
  69. if (w->symm_workspace_p)
  70. gsl_eigen_symm_free(w->symm_workspace_p);
  71. free(w);
  72. } /* gsl_eigen_gensymm_free() */
  73. /*
  74. gsl_eigen_gensymm()
  75. Solve the generalized symmetric-definite eigenvalue problem
  76. A x = \lambda B x
  77. for the eigenvalues \lambda.
  78. Inputs: A - real symmetric matrix
  79. B - real symmetric and positive definite matrix
  80. eval - where to store eigenvalues
  81. w - workspace
  82. Return: success or error
  83. */
  84. int
  85. gsl_eigen_gensymm (gsl_matrix * A, gsl_matrix * B, gsl_vector * eval,
  86. gsl_eigen_gensymm_workspace * w)
  87. {
  88. const size_t N = A->size1;
  89. /* check matrix and vector sizes */
  90. if (N != A->size2)
  91. {
  92. GSL_ERROR ("matrix must be square to compute eigenvalues", GSL_ENOTSQR);
  93. }
  94. else if ((N != B->size1) || (N != B->size2))
  95. {
  96. GSL_ERROR ("B matrix dimensions must match A", GSL_EBADLEN);
  97. }
  98. else if (eval->size != N)
  99. {
  100. GSL_ERROR ("eigenvalue vector must match matrix size", GSL_EBADLEN);
  101. }
  102. else if (w->size != N)
  103. {
  104. GSL_ERROR ("matrix size does not match workspace", GSL_EBADLEN);
  105. }
  106. else
  107. {
  108. int s;
  109. /* compute Cholesky factorization of B */
  110. s = gsl_linalg_cholesky_decomp(B);
  111. if (s != GSL_SUCCESS)
  112. return s; /* B is not positive definite */
  113. /* transform to standard symmetric eigenvalue problem */
  114. gsl_eigen_gensymm_standardize(A, B);
  115. s = gsl_eigen_symm(A, eval, w->symm_workspace_p);
  116. return s;
  117. }
  118. } /* gsl_eigen_gensymm() */
  119. /*
  120. gsl_eigen_gensymm_standardize()
  121. Reduce the generalized symmetric-definite eigenproblem to
  122. the standard symmetric eigenproblem by computing
  123. C = L^{-1} A L^{-t}
  124. where L L^t is the Cholesky decomposition of B
  125. Inputs: A - (input/output) real symmetric matrix
  126. B - real symmetric, positive definite matrix in Cholesky form
  127. Return: success
  128. Notes: A is overwritten by L^{-1} A L^{-t}
  129. */
  130. int
  131. gsl_eigen_gensymm_standardize(gsl_matrix *A, const gsl_matrix *B)
  132. {
  133. const size_t N = A->size1;
  134. size_t i;
  135. double a, b, c;
  136. for (i = 0; i < N; ++i)
  137. {
  138. /* update lower triangle of A(i:n, i:n) */
  139. a = gsl_matrix_get(A, i, i);
  140. b = gsl_matrix_get(B, i, i);
  141. a /= b * b;
  142. gsl_matrix_set(A, i, i, a);
  143. if (i < N - 1)
  144. {
  145. gsl_vector_view ai = gsl_matrix_subcolumn(A, i, i + 1, N - i - 1);
  146. gsl_matrix_view ma =
  147. gsl_matrix_submatrix(A, i + 1, i + 1, N - i - 1, N - i - 1);
  148. gsl_vector_const_view bi =
  149. gsl_matrix_const_subcolumn(B, i, i + 1, N - i - 1);
  150. gsl_matrix_const_view mb =
  151. gsl_matrix_const_submatrix(B, i + 1, i + 1, N - i - 1, N - i - 1);
  152. gsl_blas_dscal(1.0 / b, &ai.vector);
  153. c = -0.5 * a;
  154. gsl_blas_daxpy(c, &bi.vector, &ai.vector);
  155. gsl_blas_dsyr2(CblasLower, -1.0, &ai.vector, &bi.vector, &ma.matrix);
  156. gsl_blas_daxpy(c, &bi.vector, &ai.vector);
  157. gsl_blas_dtrsv(CblasLower,
  158. CblasNoTrans,
  159. CblasNonUnit,
  160. &mb.matrix,
  161. &ai.vector);
  162. }
  163. }
  164. return GSL_SUCCESS;
  165. } /* gsl_eigen_gensymm_standardize() */