gsl_linalg__hesstri.c 5.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174
  1. /* linalg/hesstri.c
  2. *
  3. * Copyright (C) 2006, 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 <math.h>
  21. #include "gsl__config.h"
  22. #include "gsl_linalg.h"
  23. #include "gsl_matrix.h"
  24. #include "gsl_vector.h"
  25. #include "gsl_blas.h"
  26. #include "gsl_linalg__givens.c"
  27. /*
  28. * This module contains routines related to the Hessenberg-Triangular
  29. * reduction of two general real matrices
  30. *
  31. * See Golub & Van Loan, "Matrix Computations", 3rd ed, sec 7.7.4
  32. */
  33. /*
  34. gsl_linalg_hesstri_decomp()
  35. Perform a reduction to generalized upper Hessenberg form.
  36. Given A and B, this function overwrites A with an upper Hessenberg
  37. matrix H = U^T A V and B with an upper triangular matrix R = U^T B V
  38. with U and V orthogonal.
  39. See Golub & Van Loan, "Matrix Computations" (3rd ed) algorithm 7.7.1
  40. Inputs: A - real square matrix
  41. B - real square matrix
  42. U - (output) if non-null, U is stored here on output
  43. V - (output) if non-null, V is stored here on output
  44. work - workspace (length n)
  45. Return: success or error
  46. */
  47. int
  48. gsl_linalg_hesstri_decomp(gsl_matrix * A, gsl_matrix * B, gsl_matrix * U,
  49. gsl_matrix * V, gsl_vector * work)
  50. {
  51. const size_t N = A->size1;
  52. if ((N != A->size2) || (N != B->size1) || (N != B->size2))
  53. {
  54. GSL_ERROR ("Hessenberg-triangular reduction requires square matrices",
  55. GSL_ENOTSQR);
  56. }
  57. else if (N != work->size)
  58. {
  59. GSL_ERROR ("length of workspace must match matrix dimension",
  60. GSL_EBADLEN);
  61. }
  62. else
  63. {
  64. double cs, sn; /* rotation parameters */
  65. size_t i, j; /* looping */
  66. gsl_vector_view xv, yv; /* temporary views */
  67. /* B -> Q^T B = R (upper triangular) */
  68. gsl_linalg_QR_decomp(B, work);
  69. /* A -> Q^T A */
  70. gsl_linalg_QR_QTmat(B, work, A);
  71. /* initialize U and V if desired */
  72. if (U)
  73. {
  74. gsl_linalg_QR_unpack(B, work, U, B);
  75. }
  76. else
  77. {
  78. /* zero out lower triangle of B */
  79. for (j = 0; j < N - 1; ++j)
  80. {
  81. for (i = j + 1; i < N; ++i)
  82. gsl_matrix_set(B, i, j, 0.0);
  83. }
  84. }
  85. if (V)
  86. gsl_matrix_set_identity(V);
  87. if (N < 3)
  88. return GSL_SUCCESS; /* nothing more to do */
  89. /* reduce A and B */
  90. for (j = 0; j < N - 2; ++j)
  91. {
  92. for (i = N - 1; i >= (j + 2); --i)
  93. {
  94. /* step 1: rotate rows i - 1, i to kill A(i,j) */
  95. /*
  96. * compute G = [ CS SN ] so that G^t [ A(i-1,j) ] = [ * ]
  97. * [-SN CS ] [ A(i, j) ] [ 0 ]
  98. */
  99. create_givens(gsl_matrix_get(A, i - 1, j),
  100. gsl_matrix_get(A, i, j),
  101. &cs,
  102. &sn);
  103. /* invert so drot() works correctly (G -> G^t) */
  104. sn = -sn;
  105. /* compute G^t A(i-1:i, j:n) */
  106. xv = gsl_matrix_subrow(A, i - 1, j, N - j);
  107. yv = gsl_matrix_subrow(A, i, j, N - j);
  108. gsl_blas_drot(&xv.vector, &yv.vector, cs, sn);
  109. /* compute G^t B(i-1:i, i-1:n) */
  110. xv = gsl_matrix_subrow(B, i - 1, i - 1, N - i + 1);
  111. yv = gsl_matrix_subrow(B, i, i - 1, N - i + 1);
  112. gsl_blas_drot(&xv.vector, &yv.vector, cs, sn);
  113. if (U)
  114. {
  115. /* accumulate U: U -> U G */
  116. xv = gsl_matrix_column(U, i - 1);
  117. yv = gsl_matrix_column(U, i);
  118. gsl_blas_drot(&xv.vector, &yv.vector, cs, sn);
  119. }
  120. /* step 2: rotate columns i, i - 1 to kill B(i, i - 1) */
  121. create_givens(-gsl_matrix_get(B, i, i),
  122. gsl_matrix_get(B, i, i - 1),
  123. &cs,
  124. &sn);
  125. /* invert so drot() works correctly (G -> G^t) */
  126. sn = -sn;
  127. /* compute B(1:i, i-1:i) G */
  128. xv = gsl_matrix_subcolumn(B, i - 1, 0, i + 1);
  129. yv = gsl_matrix_subcolumn(B, i, 0, i + 1);
  130. gsl_blas_drot(&xv.vector, &yv.vector, cs, sn);
  131. /* apply to A(1:n, i-1:i) */
  132. xv = gsl_matrix_column(A, i - 1);
  133. yv = gsl_matrix_column(A, i);
  134. gsl_blas_drot(&xv.vector, &yv.vector, cs, sn);
  135. if (V)
  136. {
  137. /* accumulate V: V -> V G */
  138. xv = gsl_matrix_column(V, i - 1);
  139. yv = gsl_matrix_column(V, i);
  140. gsl_blas_drot(&xv.vector, &yv.vector, cs, sn);
  141. }
  142. }
  143. }
  144. return GSL_SUCCESS;
  145. }
  146. } /* gsl_linalg_hesstri_decomp() */