gsl_linalg__apply_givens.c 3.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126
  1. /* linalg/apply_givens.c
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2001, 2007 Gerard Jungman, Brian Gough
  4. * Copyright (C) 2004 Joerg Wensch, modifications for LQ.
  5. *
  6. * This program is free software; you can redistribute it and/or modify
  7. * it under the terms of the GNU General Public License as published by
  8. * the Free Software Foundation; either version 3 of the License, or (at
  9. * your option) any later version.
  10. *
  11. * This program is distributed in the hope that it will be useful, but
  12. * WITHOUT ANY WARRANTY; without even the implied warranty of
  13. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
  14. * General Public License for more details.
  15. *
  16. * You should have received a copy of the GNU General Public License
  17. * along with this program; if not, write to the Free Software
  18. * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
  19. */
  20. inline static void
  21. apply_givens_qr (size_t M, size_t N, gsl_matrix * Q, gsl_matrix * R,
  22. size_t i, size_t j, double c, double s)
  23. {
  24. size_t k;
  25. /* Apply rotation to matrix Q, Q' = Q G */
  26. #if USE_BLAS
  27. {
  28. gsl_matrix_view Q0M = gsl_matrix_submatrix(Q,0,0,M,j+1);
  29. gsl_vector_view Qi = gsl_matrix_column(&Q0M.matrix,i);
  30. gsl_vector_view Qj = gsl_matrix_column(&Q0M.matrix,j);
  31. gsl_blas_drot(&Qi.vector, &Qj.vector, c, -s);
  32. }
  33. #else
  34. for (k = 0; k < M; k++)
  35. {
  36. double qki = gsl_matrix_get (Q, k, i);
  37. double qkj = gsl_matrix_get (Q, k, j);
  38. gsl_matrix_set (Q, k, i, qki * c - qkj * s);
  39. gsl_matrix_set (Q, k, j, qki * s + qkj * c);
  40. }
  41. #endif
  42. /* Apply rotation to matrix R, R' = G^T R (note: upper triangular so
  43. zero for column < row) */
  44. #if USE_BLAS
  45. {
  46. k = GSL_MIN(i,j);
  47. gsl_matrix_view R0 = gsl_matrix_submatrix(R, 0, k, j+1, N-k);
  48. gsl_vector_view Ri = gsl_matrix_row(&R0.matrix,i);
  49. gsl_vector_view Rj = gsl_matrix_row(&R0.matrix,j);
  50. gsl_blas_drot(&Ri.vector, &Rj.vector, c, -s);
  51. }
  52. #else
  53. for (k = GSL_MIN (i, j); k < N; k++)
  54. {
  55. double rik = gsl_matrix_get (R, i, k);
  56. double rjk = gsl_matrix_get (R, j, k);
  57. gsl_matrix_set (R, i, k, c * rik - s * rjk);
  58. gsl_matrix_set (R, j, k, s * rik + c * rjk);
  59. }
  60. #endif
  61. }
  62. inline static void
  63. apply_givens_lq (size_t M, size_t N, gsl_matrix * Q, gsl_matrix * L,
  64. size_t i, size_t j, double c, double s)
  65. {
  66. size_t k;
  67. /* Apply rotation to matrix Q, Q' = G Q */
  68. #if USE_BLAS
  69. {
  70. gsl_matrix_view Q0M = gsl_matrix_submatrix(Q,0,0,j+1,M);
  71. gsl_vector_view Qi = gsl_matrix_row(&Q0M.matrix,i);
  72. gsl_vector_view Qj = gsl_matrix_row(&Q0M.matrix,j);
  73. gsl_blas_drot(&Qi.vector, &Qj.vector, c, -s);
  74. }
  75. #else
  76. for (k = 0; k < M; k++)
  77. {
  78. double qik = gsl_matrix_get (Q, i, k);
  79. double qjk = gsl_matrix_get (Q, j, k);
  80. gsl_matrix_set (Q, i, k, qik * c - qjk * s);
  81. gsl_matrix_set (Q, j, k, qik * s + qjk * c);
  82. }
  83. #endif
  84. /* Apply rotation to matrix L, L' = L G^T (note: lower triangular so
  85. zero for column > row) */
  86. #if USE_BLAS
  87. {
  88. k = GSL_MIN(i,j);
  89. gsl_matrix_view L0 = gsl_matrix_submatrix(L, k, 0, N-k, j+1);
  90. gsl_vector_view Li = gsl_matrix_column(&L0.matrix,i);
  91. gsl_vector_view Lj = gsl_matrix_column(&L0.matrix,j);
  92. gsl_blas_drot(&Li.vector, &Lj.vector, c, -s);
  93. }
  94. #else
  95. for (k = GSL_MIN (i, j); k < N; k++)
  96. {
  97. double lki = gsl_matrix_get (L, k, i);
  98. double lkj = gsl_matrix_get (L, k, j);
  99. gsl_matrix_set (L, k, i, c * lki - s * lkj);
  100. gsl_matrix_set (L, k, j, s * lki + c * lkj);
  101. }
  102. #endif
  103. }
  104. inline static void
  105. apply_givens_vec (gsl_vector * v, size_t i, size_t j, double c, double s)
  106. {
  107. /* Apply rotation to vector v' = G^T v */
  108. double vi = gsl_vector_get (v, i);
  109. double vj = gsl_vector_get (v, j);
  110. gsl_vector_set (v, i, c * vi - s * vj);
  111. gsl_vector_set (v, j, s * vi + c * vj);
  112. }