gsl_linalg__hh.c 4.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180
  1. /* linalg/hh.c
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000, 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. /* Author: G. Jungman */
  20. #include "gsl__config.h"
  21. #include <stdlib.h>
  22. #include "gsl_math.h"
  23. #include "gsl_vector.h"
  24. #include "gsl_matrix.h"
  25. #include "gsl_linalg.h"
  26. #define REAL double
  27. /* [Engeln-Mullges + Uhlig, Alg. 4.42]
  28. */
  29. int
  30. gsl_linalg_HH_solve (gsl_matrix * A, const gsl_vector * b, gsl_vector * x)
  31. {
  32. if (A->size1 > A->size2)
  33. {
  34. /* System is underdetermined. */
  35. GSL_ERROR ("System is underdetermined", GSL_EINVAL);
  36. }
  37. else if (A->size2 != x->size)
  38. {
  39. GSL_ERROR ("matrix and vector sizes must be equal", GSL_EBADLEN);
  40. }
  41. else
  42. {
  43. int status ;
  44. gsl_vector_memcpy (x, b);
  45. status = gsl_linalg_HH_svx (A, x);
  46. return status ;
  47. }
  48. }
  49. int
  50. gsl_linalg_HH_svx (gsl_matrix * A, gsl_vector * x)
  51. {
  52. if (A->size1 > A->size2)
  53. {
  54. /* System is underdetermined. */
  55. GSL_ERROR ("System is underdetermined", GSL_EINVAL);
  56. }
  57. else if (A->size2 != x->size)
  58. {
  59. GSL_ERROR ("matrix and vector sizes must be equal", GSL_EBADLEN);
  60. }
  61. else
  62. {
  63. const size_t N = A->size1;
  64. const size_t M = A->size2;
  65. size_t i, j, k;
  66. REAL *d = (REAL *) malloc (N * sizeof (REAL));
  67. if (d == 0)
  68. {
  69. GSL_ERROR ("could not allocate memory for workspace", GSL_ENOMEM);
  70. }
  71. /* Perform Householder transformation. */
  72. for (i = 0; i < N; i++)
  73. {
  74. const REAL aii = gsl_matrix_get (A, i, i);
  75. REAL alpha;
  76. REAL f;
  77. REAL ak;
  78. REAL max_norm = 0.0;
  79. REAL r = 0.0;
  80. for (k = i; k < M; k++)
  81. {
  82. REAL aki = gsl_matrix_get (A, k, i);
  83. r += aki * aki;
  84. }
  85. if (r == 0.0)
  86. {
  87. /* Rank of matrix is less than size1. */
  88. free (d);
  89. GSL_ERROR ("matrix is rank deficient", GSL_ESING);
  90. }
  91. alpha = sqrt (r) * GSL_SIGN (aii);
  92. ak = 1.0 / (r + alpha * aii);
  93. gsl_matrix_set (A, i, i, aii + alpha);
  94. d[i] = -alpha;
  95. for (k = i + 1; k < N; k++)
  96. {
  97. REAL norm = 0.0;
  98. f = 0.0;
  99. for (j = i; j < M; j++)
  100. {
  101. REAL ajk = gsl_matrix_get (A, j, k);
  102. REAL aji = gsl_matrix_get (A, j, i);
  103. norm += ajk * ajk;
  104. f += ajk * aji;
  105. }
  106. max_norm = GSL_MAX (max_norm, norm);
  107. f *= ak;
  108. for (j = i; j < M; j++)
  109. {
  110. REAL ajk = gsl_matrix_get (A, j, k);
  111. REAL aji = gsl_matrix_get (A, j, i);
  112. gsl_matrix_set (A, j, k, ajk - f * aji);
  113. }
  114. }
  115. if (fabs (alpha) < 2.0 * GSL_DBL_EPSILON * sqrt (max_norm))
  116. {
  117. /* Apparent singularity. */
  118. free (d);
  119. GSL_ERROR("apparent singularity detected", GSL_ESING);
  120. }
  121. /* Perform update of RHS. */
  122. f = 0.0;
  123. for (j = i; j < M; j++)
  124. {
  125. f += gsl_vector_get (x, j) * gsl_matrix_get (A, j, i);
  126. }
  127. f *= ak;
  128. for (j = i; j < M; j++)
  129. {
  130. REAL xj = gsl_vector_get (x, j);
  131. REAL aji = gsl_matrix_get (A, j, i);
  132. gsl_vector_set (x, j, xj - f * aji);
  133. }
  134. }
  135. /* Perform back-substitution. */
  136. for (i = N; i > 0 && i--;)
  137. {
  138. REAL xi = gsl_vector_get (x, i);
  139. REAL sum = 0.0;
  140. for (k = i + 1; k < N; k++)
  141. {
  142. sum += gsl_matrix_get (A, i, k) * gsl_vector_get (x, k);
  143. }
  144. gsl_vector_set (x, i, (xi - sum) / d[i]);
  145. }
  146. free (d);
  147. return GSL_SUCCESS;
  148. }
  149. }