gsl_linalg__householdercomplex.c 7.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261
  1. /* linalg/householdercomplex.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. /* Computes a householder transformation matrix H such that
  20. *
  21. * H' v = -/+ |v| e_1
  22. *
  23. * where e_1 is the first unit vector. On exit the matrix H can be
  24. * computed from the return values (tau, v)
  25. *
  26. * H = I - tau * w * w'
  27. *
  28. * where w = (1, v(2), ..., v(N)). The nonzero element of the result
  29. * vector -/+|v| e_1 is stored in v(1).
  30. *
  31. * Note that the matrix H' in the householder transformation is the
  32. * hermitian conjugate of H. To compute H'v, pass the conjugate of
  33. * tau as the first argument to gsl_linalg_householder_hm() rather
  34. * than tau itself. See the LAPACK function CLARFG for details of this
  35. * convention. */
  36. #include "gsl__config.h"
  37. #include <stdlib.h>
  38. #include "gsl_math.h"
  39. #include "gsl_vector.h"
  40. #include "gsl_matrix.h"
  41. #include "gsl_blas.h"
  42. #include "gsl_complex_math.h"
  43. #include "gsl_linalg.h"
  44. gsl_complex
  45. gsl_linalg_complex_householder_transform (gsl_vector_complex * v)
  46. {
  47. /* replace v[0:n-1] with a householder vector (v[0:n-1]) and
  48. coefficient tau that annihilate v[1:n-1] */
  49. const size_t n = v->size ;
  50. if (n == 1)
  51. {
  52. gsl_complex alpha = gsl_vector_complex_get (v, 0) ;
  53. double absa = gsl_complex_abs (alpha);
  54. double beta_r = - (GSL_REAL(alpha) >= 0 ? +1 : -1) * absa ;
  55. gsl_complex tau;
  56. if (beta_r == 0.0)
  57. {
  58. GSL_REAL(tau) = 0.0;
  59. GSL_IMAG(tau) = 0.0;
  60. }
  61. else
  62. {
  63. GSL_REAL(tau) = (beta_r - GSL_REAL(alpha)) / beta_r ;
  64. GSL_IMAG(tau) = - GSL_IMAG(alpha) / beta_r ;
  65. {
  66. gsl_complex beta = gsl_complex_rect (beta_r, 0.0);
  67. gsl_vector_complex_set (v, 0, beta) ;
  68. }
  69. }
  70. return tau;
  71. }
  72. else
  73. {
  74. gsl_complex tau ;
  75. double beta_r;
  76. gsl_vector_complex_view x = gsl_vector_complex_subvector (v, 1, n - 1) ;
  77. gsl_complex alpha = gsl_vector_complex_get (v, 0) ;
  78. double absa = gsl_complex_abs (alpha);
  79. double xnorm = gsl_blas_dznrm2 (&x.vector);
  80. if (xnorm == 0 && GSL_IMAG(alpha) == 0)
  81. {
  82. gsl_complex zero = gsl_complex_rect(0.0, 0.0);
  83. return zero; /* tau = 0 */
  84. }
  85. beta_r = - (GSL_REAL(alpha) >= 0 ? +1 : -1) * hypot(absa, xnorm) ;
  86. GSL_REAL(tau) = (beta_r - GSL_REAL(alpha)) / beta_r ;
  87. GSL_IMAG(tau) = - GSL_IMAG(alpha) / beta_r ;
  88. {
  89. gsl_complex amb = gsl_complex_sub_real(alpha, beta_r);
  90. gsl_complex s = gsl_complex_inverse(amb);
  91. gsl_blas_zscal (s, &x.vector);
  92. }
  93. {
  94. gsl_complex beta = gsl_complex_rect (beta_r, 0.0);
  95. gsl_vector_complex_set (v, 0, beta) ;
  96. }
  97. return tau;
  98. }
  99. }
  100. int
  101. gsl_linalg_complex_householder_hm (gsl_complex tau, const gsl_vector_complex * v, gsl_matrix_complex * A)
  102. {
  103. /* applies a householder transformation v,tau to matrix m */
  104. size_t i, j;
  105. if (GSL_REAL(tau) == 0.0 && GSL_IMAG(tau) == 0.0)
  106. {
  107. return GSL_SUCCESS;
  108. }
  109. /* w = (v' A)^T */
  110. for (j = 0; j < A->size2; j++)
  111. {
  112. gsl_complex tauwj;
  113. gsl_complex wj = gsl_matrix_complex_get(A,0,j);
  114. for (i = 1; i < A->size1; i++) /* note, computed for v(0) = 1 above */
  115. {
  116. gsl_complex Aij = gsl_matrix_complex_get(A,i,j);
  117. gsl_complex vi = gsl_vector_complex_get(v,i);
  118. gsl_complex Av = gsl_complex_mul (Aij, gsl_complex_conjugate(vi));
  119. wj = gsl_complex_add (wj, Av);
  120. }
  121. tauwj = gsl_complex_mul (tau, wj);
  122. /* A = A - v w^T */
  123. {
  124. gsl_complex A0j = gsl_matrix_complex_get (A, 0, j);
  125. gsl_complex Atw = gsl_complex_sub (A0j, tauwj);
  126. /* store A0j - tau * wj */
  127. gsl_matrix_complex_set (A, 0, j, Atw);
  128. }
  129. for (i = 1; i < A->size1; i++)
  130. {
  131. gsl_complex vi = gsl_vector_complex_get (v, i);
  132. gsl_complex tauvw = gsl_complex_mul(vi, tauwj);
  133. gsl_complex Aij = gsl_matrix_complex_get (A, i, j);
  134. gsl_complex Atwv = gsl_complex_sub (Aij, tauvw);
  135. /* store Aij - tau * vi * wj */
  136. gsl_matrix_complex_set (A, i, j, Atwv);
  137. }
  138. }
  139. return GSL_SUCCESS;
  140. }
  141. int
  142. gsl_linalg_complex_householder_mh (gsl_complex tau, const gsl_vector_complex * v, gsl_matrix_complex * A)
  143. {
  144. /* applies a householder transformation v,tau to matrix m on the right */
  145. size_t i, j;
  146. if (GSL_REAL(tau) == 0.0 && GSL_IMAG(tau) == 0.0)
  147. {
  148. return GSL_SUCCESS;
  149. }
  150. /* A -> A - A*tau*v*v^h */
  151. for (i = 0; i < A->size1; i++)
  152. {
  153. gsl_complex tauwi;
  154. gsl_complex Ai0 = gsl_matrix_complex_get (A, i, 0);
  155. gsl_complex wi = Ai0;
  156. /* compute w = A v */
  157. for (j = 1; j < A->size2; j++) /* note, computed for v(0) = 1 above */
  158. {
  159. gsl_complex Aij = gsl_matrix_complex_get(A, i, j);
  160. gsl_complex vj = gsl_vector_complex_get(v, j);
  161. gsl_complex Av = gsl_complex_mul (Aij, vj);
  162. wi = gsl_complex_add (wi, Av);
  163. }
  164. tauwi = gsl_complex_mul (tau, wi);
  165. /* A = A - w v^H */
  166. {
  167. gsl_complex Atw = gsl_complex_sub (Ai0, tauwi);
  168. /* store Ai0 - tau * wi */
  169. gsl_matrix_complex_set (A, i, 0, Atw);
  170. }
  171. for (j = 1; j < A->size2; j++)
  172. {
  173. gsl_complex vj = gsl_vector_complex_get (v, j);
  174. gsl_complex tauwv = gsl_complex_mul(gsl_complex_conjugate(vj), tauwi);
  175. gsl_complex Aij = gsl_matrix_complex_get (A, i, j);
  176. gsl_complex Atwv = gsl_complex_sub (Aij, tauwv);
  177. /* store Aij - tau * wi * conj(vj) */
  178. gsl_matrix_complex_set (A, i, j, Atwv);
  179. }
  180. }
  181. return GSL_SUCCESS;
  182. }
  183. int
  184. gsl_linalg_complex_householder_hv (gsl_complex tau, const gsl_vector_complex * v, gsl_vector_complex * w)
  185. {
  186. const size_t N = v->size;
  187. if (GSL_REAL(tau) == 0.0 && GSL_IMAG(tau) == 0.0)
  188. return GSL_SUCCESS;
  189. {
  190. /* compute z = v'w */
  191. gsl_complex z0 = gsl_vector_complex_get(w,0);
  192. gsl_complex z1, z;
  193. gsl_complex tz, ntz;
  194. gsl_vector_complex_const_view v1 = gsl_vector_complex_const_subvector(v, 1, N-1);
  195. gsl_vector_complex_view w1 = gsl_vector_complex_subvector(w, 1, N-1);
  196. gsl_blas_zdotc(&v1.vector, &w1.vector, &z1);
  197. z = gsl_complex_add (z0, z1);
  198. tz = gsl_complex_mul(tau, z);
  199. ntz = gsl_complex_negative (tz);
  200. /* compute w = w - tau * (v'w) * v */
  201. {
  202. gsl_complex w0 = gsl_vector_complex_get(w, 0);
  203. gsl_complex w0ntz = gsl_complex_add (w0, ntz);
  204. gsl_vector_complex_set (w, 0, w0ntz);
  205. }
  206. gsl_blas_zaxpy(ntz, &v1.vector, &w1.vector);
  207. }
  208. return GSL_SUCCESS;
  209. }