gsl_linalg__householder.c 7.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327
  1. /* linalg/householder.c
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2004, 2007 Gerard Jungman, 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. #include "gsl__config.h"
  20. #include <stdlib.h>
  21. #include "gsl_math.h"
  22. #include "gsl_vector.h"
  23. #include "gsl_matrix.h"
  24. #include "gsl_blas.h"
  25. #include "gsl_linalg.h"
  26. double
  27. gsl_linalg_householder_transform (gsl_vector * v)
  28. {
  29. /* replace v[0:n-1] with a householder vector (v[0:n-1]) and
  30. coefficient tau that annihilate v[1:n-1] */
  31. const size_t n = v->size ;
  32. if (n == 1)
  33. {
  34. return 0.0; /* tau = 0 */
  35. }
  36. else
  37. {
  38. double alpha, beta, tau ;
  39. gsl_vector_view x = gsl_vector_subvector (v, 1, n - 1) ;
  40. double xnorm = gsl_blas_dnrm2 (&x.vector);
  41. if (xnorm == 0)
  42. {
  43. return 0.0; /* tau = 0 */
  44. }
  45. alpha = gsl_vector_get (v, 0) ;
  46. beta = - (alpha >= 0.0 ? +1.0 : -1.0) * hypot(alpha, xnorm) ;
  47. tau = (beta - alpha) / beta ;
  48. gsl_blas_dscal (1.0 / (alpha - beta), &x.vector);
  49. gsl_vector_set (v, 0, beta) ;
  50. return tau;
  51. }
  52. }
  53. int
  54. gsl_linalg_householder_hm (double tau, const gsl_vector * v, gsl_matrix * A)
  55. {
  56. /* applies a householder transformation v,tau to matrix m */
  57. if (tau == 0.0)
  58. {
  59. return GSL_SUCCESS;
  60. }
  61. #ifdef USE_BLAS
  62. {
  63. gsl_vector_const_view v1 = gsl_vector_const_subvector (v, 1, v->size - 1);
  64. gsl_matrix_view A1 = gsl_matrix_submatrix (A, 1, 0, A->size1 - 1, A->size2);
  65. size_t j;
  66. for (j = 0; j < A->size2; j++)
  67. {
  68. double wj = 0.0;
  69. gsl_vector_view A1j = gsl_matrix_column(&A1.matrix, j);
  70. gsl_blas_ddot (&A1j.vector, &v1.vector, &wj);
  71. wj += gsl_matrix_get(A,0,j);
  72. {
  73. double A0j = gsl_matrix_get (A, 0, j);
  74. gsl_matrix_set (A, 0, j, A0j - tau * wj);
  75. }
  76. gsl_blas_daxpy (-tau * wj, &v1.vector, &A1j.vector);
  77. }
  78. }
  79. #else
  80. {
  81. size_t i, j;
  82. for (j = 0; j < A->size2; j++)
  83. {
  84. /* Compute wj = Akj vk */
  85. double wj = gsl_matrix_get(A,0,j);
  86. for (i = 1; i < A->size1; i++) /* note, computed for v(0) = 1 above */
  87. {
  88. wj += gsl_matrix_get(A,i,j) * gsl_vector_get(v,i);
  89. }
  90. /* Aij = Aij - tau vi wj */
  91. /* i = 0 */
  92. {
  93. double A0j = gsl_matrix_get (A, 0, j);
  94. gsl_matrix_set (A, 0, j, A0j - tau * wj);
  95. }
  96. /* i = 1 .. M-1 */
  97. for (i = 1; i < A->size1; i++)
  98. {
  99. double Aij = gsl_matrix_get (A, i, j);
  100. double vi = gsl_vector_get (v, i);
  101. gsl_matrix_set (A, i, j, Aij - tau * vi * wj);
  102. }
  103. }
  104. }
  105. #endif
  106. return GSL_SUCCESS;
  107. }
  108. int
  109. gsl_linalg_householder_mh (double tau, const gsl_vector * v, gsl_matrix * A)
  110. {
  111. /* applies a householder transformation v,tau to matrix m from the
  112. right hand side in order to zero out rows */
  113. if (tau == 0)
  114. return GSL_SUCCESS;
  115. /* A = A - tau w v' */
  116. #ifdef USE_BLAS
  117. {
  118. gsl_vector_const_view v1 = gsl_vector_const_subvector (v, 1, v->size - 1);
  119. gsl_matrix_view A1 = gsl_matrix_submatrix (A, 0, 1, A->size1, A->size2-1);
  120. size_t i;
  121. for (i = 0; i < A->size1; i++)
  122. {
  123. double wi = 0.0;
  124. gsl_vector_view A1i = gsl_matrix_row(&A1.matrix, i);
  125. gsl_blas_ddot (&A1i.vector, &v1.vector, &wi);
  126. wi += gsl_matrix_get(A,i,0);
  127. {
  128. double Ai0 = gsl_matrix_get (A, i, 0);
  129. gsl_matrix_set (A, i, 0, Ai0 - tau * wi);
  130. }
  131. gsl_blas_daxpy(-tau * wi, &v1.vector, &A1i.vector);
  132. }
  133. }
  134. #else
  135. {
  136. size_t i, j;
  137. for (i = 0; i < A->size1; i++)
  138. {
  139. double wi = gsl_matrix_get(A,i,0);
  140. for (j = 1; j < A->size2; j++) /* note, computed for v(0) = 1 above */
  141. {
  142. wi += gsl_matrix_get(A,i,j) * gsl_vector_get(v,j);
  143. }
  144. /* j = 0 */
  145. {
  146. double Ai0 = gsl_matrix_get (A, i, 0);
  147. gsl_matrix_set (A, i, 0, Ai0 - tau * wi);
  148. }
  149. /* j = 1 .. N-1 */
  150. for (j = 1; j < A->size2; j++)
  151. {
  152. double vj = gsl_vector_get (v, j);
  153. double Aij = gsl_matrix_get (A, i, j);
  154. gsl_matrix_set (A, i, j, Aij - tau * wi * vj);
  155. }
  156. }
  157. }
  158. #endif
  159. return GSL_SUCCESS;
  160. }
  161. int
  162. gsl_linalg_householder_hv (double tau, const gsl_vector * v, gsl_vector * w)
  163. {
  164. /* applies a householder transformation v to vector w */
  165. const size_t N = v->size;
  166. if (tau == 0)
  167. return GSL_SUCCESS ;
  168. {
  169. /* compute d = v'w */
  170. double d0 = gsl_vector_get(w,0);
  171. double d1, d;
  172. gsl_vector_const_view v1 = gsl_vector_const_subvector(v, 1, N-1);
  173. gsl_vector_view w1 = gsl_vector_subvector(w, 1, N-1);
  174. gsl_blas_ddot (&v1.vector, &w1.vector, &d1);
  175. d = d0 + d1;
  176. /* compute w = w - tau (v) (v'w) */
  177. {
  178. double w0 = gsl_vector_get (w,0);
  179. gsl_vector_set (w, 0, w0 - tau * d);
  180. }
  181. gsl_blas_daxpy (-tau * d, &v1.vector, &w1.vector);
  182. }
  183. return GSL_SUCCESS;
  184. }
  185. int
  186. gsl_linalg_householder_hm1 (double tau, gsl_matrix * A)
  187. {
  188. /* applies a householder transformation v,tau to a matrix being
  189. build up from the identity matrix, using the first column of A as
  190. a householder vector */
  191. if (tau == 0)
  192. {
  193. size_t i,j;
  194. gsl_matrix_set (A, 0, 0, 1.0);
  195. for (j = 1; j < A->size2; j++)
  196. {
  197. gsl_matrix_set (A, 0, j, 0.0);
  198. }
  199. for (i = 1; i < A->size1; i++)
  200. {
  201. gsl_matrix_set (A, i, 0, 0.0);
  202. }
  203. return GSL_SUCCESS;
  204. }
  205. /* w = A' v */
  206. #ifdef USE_BLAS
  207. {
  208. gsl_matrix_view A1 = gsl_matrix_submatrix (A, 1, 0, A->size1 - 1, A->size2);
  209. gsl_vector_view v1 = gsl_matrix_column (&A1.matrix, 0);
  210. size_t j;
  211. for (j = 1; j < A->size2; j++)
  212. {
  213. double wj = 0.0; /* A0j * v0 */
  214. gsl_vector_view A1j = gsl_matrix_column(&A1.matrix, j);
  215. gsl_blas_ddot (&A1j.vector, &v1.vector, &wj);
  216. /* A = A - tau v w' */
  217. gsl_matrix_set (A, 0, j, - tau * wj);
  218. gsl_blas_daxpy(-tau*wj, &v1.vector, &A1j.vector);
  219. }
  220. gsl_blas_dscal(-tau, &v1.vector);
  221. gsl_matrix_set (A, 0, 0, 1.0 - tau);
  222. }
  223. #else
  224. {
  225. size_t i, j;
  226. for (j = 1; j < A->size2; j++)
  227. {
  228. double wj = 0.0; /* A0j * v0 */
  229. for (i = 1; i < A->size1; i++)
  230. {
  231. double vi = gsl_matrix_get(A, i, 0);
  232. wj += gsl_matrix_get(A,i,j) * vi;
  233. }
  234. /* A = A - tau v w' */
  235. gsl_matrix_set (A, 0, j, - tau * wj);
  236. for (i = 1; i < A->size1; i++)
  237. {
  238. double vi = gsl_matrix_get (A, i, 0);
  239. double Aij = gsl_matrix_get (A, i, j);
  240. gsl_matrix_set (A, i, j, Aij - tau * vi * wj);
  241. }
  242. }
  243. for (i = 1; i < A->size1; i++)
  244. {
  245. double vi = gsl_matrix_get(A, i, 0);
  246. gsl_matrix_set(A, i, 0, -tau * vi);
  247. }
  248. gsl_matrix_set (A, 0, 0, 1.0 - tau);
  249. }
  250. #endif
  251. return GSL_SUCCESS;
  252. }