gsl_linalg__symmtd.c 6.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233
  1. /* linalg/sytd.c
  2. *
  3. * Copyright (C) 2001, 2007 Brian Gough
  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. /* Factorise a symmetric matrix A into
  20. *
  21. * A = Q T Q'
  22. *
  23. * where Q is orthogonal and T is symmetric tridiagonal. Only the
  24. * diagonal and lower triangular part of A is referenced and modified.
  25. *
  26. * On exit, T is stored in the diagonal and first subdiagonal of
  27. * A. Since T is symmetric the upper diagonal is not stored.
  28. *
  29. * Q is stored as a packed set of Householder transformations in the
  30. * lower triangular part of the input matrix below the first subdiagonal.
  31. *
  32. * The full matrix for Q can be obtained as the product
  33. *
  34. * Q = Q_1 Q_2 ... Q_(N-2)
  35. *
  36. * where
  37. *
  38. * Q_i = (I - tau_i * v_i * v_i')
  39. *
  40. * and where v_i is a Householder vector
  41. *
  42. * v_i = [0, ... , 0, 1, A(i+1,i), A(i+2,i), ... , A(N,i)]
  43. *
  44. * This storage scheme is the same as in LAPACK. See LAPACK's
  45. * ssytd2.f for details.
  46. *
  47. * See Golub & Van Loan, "Matrix Computations" (3rd ed), Section 8.3
  48. *
  49. * Note: this description uses 1-based indices. The code below uses
  50. * 0-based indices
  51. */
  52. #include "gsl__config.h"
  53. #include <stdlib.h>
  54. #include "gsl_math.h"
  55. #include "gsl_vector.h"
  56. #include "gsl_matrix.h"
  57. #include "gsl_blas.h"
  58. #include "gsl_linalg.h"
  59. int
  60. gsl_linalg_symmtd_decomp (gsl_matrix * A, gsl_vector * tau)
  61. {
  62. if (A->size1 != A->size2)
  63. {
  64. GSL_ERROR ("symmetric tridiagonal decomposition requires square matrix",
  65. GSL_ENOTSQR);
  66. }
  67. else if (tau->size + 1 != A->size1)
  68. {
  69. GSL_ERROR ("size of tau must be (matrix size - 1)", GSL_EBADLEN);
  70. }
  71. else
  72. {
  73. const size_t N = A->size1;
  74. size_t i;
  75. for (i = 0 ; i < N - 2; i++)
  76. {
  77. gsl_vector_view c = gsl_matrix_column (A, i);
  78. gsl_vector_view v = gsl_vector_subvector (&c.vector, i + 1, N - (i + 1));
  79. double tau_i = gsl_linalg_householder_transform (&v.vector);
  80. /* Apply the transformation H^T A H to the remaining columns */
  81. if (tau_i != 0.0)
  82. {
  83. gsl_matrix_view m = gsl_matrix_submatrix (A, i + 1, i + 1,
  84. N - (i+1), N - (i+1));
  85. double ei = gsl_vector_get(&v.vector, 0);
  86. gsl_vector_view x = gsl_vector_subvector (tau, i, N-(i+1));
  87. gsl_vector_set (&v.vector, 0, 1.0);
  88. /* x = tau * A * v */
  89. gsl_blas_dsymv (CblasLower, tau_i, &m.matrix, &v.vector, 0.0, &x.vector);
  90. /* w = x - (1/2) tau * (x' * v) * v */
  91. {
  92. double xv, alpha;
  93. gsl_blas_ddot(&x.vector, &v.vector, &xv);
  94. alpha = - (tau_i / 2.0) * xv;
  95. gsl_blas_daxpy(alpha, &v.vector, &x.vector);
  96. }
  97. /* apply the transformation A = A - v w' - w v' */
  98. gsl_blas_dsyr2(CblasLower, -1.0, &v.vector, &x.vector, &m.matrix);
  99. gsl_vector_set (&v.vector, 0, ei);
  100. }
  101. gsl_vector_set (tau, i, tau_i);
  102. }
  103. return GSL_SUCCESS;
  104. }
  105. }
  106. /* Form the orthogonal matrix Q from the packed QR matrix */
  107. int
  108. gsl_linalg_symmtd_unpack (const gsl_matrix * A,
  109. const gsl_vector * tau,
  110. gsl_matrix * Q,
  111. gsl_vector * diag,
  112. gsl_vector * sdiag)
  113. {
  114. if (A->size1 != A->size2)
  115. {
  116. GSL_ERROR ("matrix A must be square", GSL_ENOTSQR);
  117. }
  118. else if (tau->size + 1 != A->size1)
  119. {
  120. GSL_ERROR ("size of tau must be (matrix size - 1)", GSL_EBADLEN);
  121. }
  122. else if (Q->size1 != A->size1 || Q->size2 != A->size1)
  123. {
  124. GSL_ERROR ("size of Q must match size of A", GSL_EBADLEN);
  125. }
  126. else if (diag->size != A->size1)
  127. {
  128. GSL_ERROR ("size of diagonal must match size of A", GSL_EBADLEN);
  129. }
  130. else if (sdiag->size + 1 != A->size1)
  131. {
  132. GSL_ERROR ("size of subdiagonal must be (matrix size - 1)", GSL_EBADLEN);
  133. }
  134. else
  135. {
  136. const size_t N = A->size1;
  137. size_t i;
  138. /* Initialize Q to the identity */
  139. gsl_matrix_set_identity (Q);
  140. for (i = N - 2; i > 0 && i--;)
  141. {
  142. gsl_vector_const_view c = gsl_matrix_const_column (A, i);
  143. gsl_vector_const_view h = gsl_vector_const_subvector (&c.vector, i + 1, N - (i+1));
  144. double ti = gsl_vector_get (tau, i);
  145. gsl_matrix_view m = gsl_matrix_submatrix (Q, i + 1, i + 1, N-(i+1), N-(i+1));
  146. gsl_linalg_householder_hm (ti, &h.vector, &m.matrix);
  147. }
  148. /* Copy diagonal into diag */
  149. for (i = 0; i < N; i++)
  150. {
  151. double Aii = gsl_matrix_get (A, i, i);
  152. gsl_vector_set (diag, i, Aii);
  153. }
  154. /* Copy subdiagonal into sd */
  155. for (i = 0; i < N - 1; i++)
  156. {
  157. double Aji = gsl_matrix_get (A, i+1, i);
  158. gsl_vector_set (sdiag, i, Aji);
  159. }
  160. return GSL_SUCCESS;
  161. }
  162. }
  163. int
  164. gsl_linalg_symmtd_unpack_T (const gsl_matrix * A,
  165. gsl_vector * diag,
  166. gsl_vector * sdiag)
  167. {
  168. if (A->size1 != A->size2)
  169. {
  170. GSL_ERROR ("matrix A must be square", GSL_ENOTSQR);
  171. }
  172. else if (diag->size != A->size1)
  173. {
  174. GSL_ERROR ("size of diagonal must match size of A", GSL_EBADLEN);
  175. }
  176. else if (sdiag->size + 1 != A->size1)
  177. {
  178. GSL_ERROR ("size of subdiagonal must be (matrix size - 1)", GSL_EBADLEN);
  179. }
  180. else
  181. {
  182. const size_t N = A->size1;
  183. size_t i;
  184. /* Copy diagonal into diag */
  185. for (i = 0; i < N; i++)
  186. {
  187. double Aii = gsl_matrix_get (A, i, i);
  188. gsl_vector_set (diag, i, Aii);
  189. }
  190. /* Copy subdiagonal into sdiag */
  191. for (i = 0; i < N - 1; i++)
  192. {
  193. double Aij = gsl_matrix_get (A, i+1, i);
  194. gsl_vector_set (sdiag, i, Aij);
  195. }
  196. return GSL_SUCCESS;
  197. }
  198. }