gsl_linalg__balancemat.c 4.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187
  1. /* linalg/balance.c
  2. *
  3. * Copyright (C) 2006 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. /* Balance a general matrix by scaling the rows and columns, so the
  20. * new row and column norms are the same order of magnitude.
  21. *
  22. * B = D^-1 A D
  23. *
  24. * where D is a diagonal matrix
  25. *
  26. * This is necessary for the unsymmetric eigenvalue problem since the
  27. * calculation can become numerically unstable for unbalanced
  28. * matrices.
  29. *
  30. * See Golub & Van Loan, "Matrix Computations" (3rd ed), Section 7.5.7
  31. * and Wilkinson & Reinsch, "Handbook for Automatic Computation", II/11 p320.
  32. */
  33. #include "gsl__config.h"
  34. #include <stdlib.h>
  35. #include "gsl_math.h"
  36. #include "gsl_vector.h"
  37. #include "gsl_matrix.h"
  38. #include "gsl_blas.h"
  39. #include "gsl_linalg.h"
  40. #define FLOAT_RADIX 2.0
  41. #define FLOAT_RADIX_SQ (FLOAT_RADIX * FLOAT_RADIX)
  42. int
  43. gsl_linalg_balance_matrix(gsl_matrix * A, gsl_vector * D)
  44. {
  45. const size_t N = A->size1;
  46. if (N != D->size)
  47. {
  48. GSL_ERROR ("vector must match matrix size", GSL_EBADLEN);
  49. }
  50. else
  51. {
  52. double row_norm,
  53. col_norm;
  54. int not_converged;
  55. gsl_vector_view v;
  56. /* initialize D to the identity matrix */
  57. gsl_vector_set_all(D, 1.0);
  58. not_converged = 1;
  59. while (not_converged)
  60. {
  61. size_t i, j;
  62. double g, f, s;
  63. not_converged = 0;
  64. for (i = 0; i < N; ++i)
  65. {
  66. row_norm = 0.0;
  67. col_norm = 0.0;
  68. for (j = 0; j < N; ++j)
  69. {
  70. if (j != i)
  71. {
  72. col_norm += fabs(gsl_matrix_get(A, j, i));
  73. row_norm += fabs(gsl_matrix_get(A, i, j));
  74. }
  75. }
  76. if ((col_norm == 0.0) || (row_norm == 0.0))
  77. {
  78. continue;
  79. }
  80. g = row_norm / FLOAT_RADIX;
  81. f = 1.0;
  82. s = col_norm + row_norm;
  83. /*
  84. * find the integer power of the machine radix which
  85. * comes closest to balancing the matrix
  86. */
  87. while (col_norm < g)
  88. {
  89. f *= FLOAT_RADIX;
  90. col_norm *= FLOAT_RADIX_SQ;
  91. }
  92. g = row_norm * FLOAT_RADIX;
  93. while (col_norm > g)
  94. {
  95. f /= FLOAT_RADIX;
  96. col_norm /= FLOAT_RADIX_SQ;
  97. }
  98. if ((row_norm + col_norm) < 0.95 * s * f)
  99. {
  100. not_converged = 1;
  101. g = 1.0 / f;
  102. /*
  103. * apply similarity transformation D, where
  104. * D_{ij} = f_i * delta_{ij}
  105. */
  106. /* multiply by D^{-1} on the left */
  107. v = gsl_matrix_row(A, i);
  108. gsl_blas_dscal(g, &v.vector);
  109. /* multiply by D on the right */
  110. v = gsl_matrix_column(A, i);
  111. gsl_blas_dscal(f, &v.vector);
  112. /* keep track of transformation */
  113. gsl_vector_set(D, i, gsl_vector_get(D, i) * f);
  114. }
  115. }
  116. }
  117. return GSL_SUCCESS;
  118. }
  119. } /* gsl_linalg_balance_matrix() */
  120. /*
  121. gsl_linalg_balance_accum()
  122. Accumulate a balancing transformation into a matrix.
  123. This is used during the computation of Schur vectors since the
  124. Schur vectors computed are the vectors for the balanced matrix.
  125. We must at some point accumulate the balancing transformation into
  126. the Schur vector matrix to get the vectors for the original matrix.
  127. A -> D A
  128. where D is the diagonal matrix
  129. Inputs: A - matrix to transform
  130. D - vector containing diagonal elements of D
  131. */
  132. int
  133. gsl_linalg_balance_accum(gsl_matrix *A, gsl_vector *D)
  134. {
  135. const size_t N = A->size1;
  136. if (N != D->size)
  137. {
  138. GSL_ERROR ("vector must match matrix size", GSL_EBADLEN);
  139. }
  140. else
  141. {
  142. size_t i;
  143. double s;
  144. gsl_vector_view r;
  145. for (i = 0; i < N; ++i)
  146. {
  147. s = gsl_vector_get(D, i);
  148. r = gsl_matrix_row(A, i);
  149. gsl_blas_dscal(s, &r.vector);
  150. }
  151. return GSL_SUCCESS;
  152. }
  153. } /* gsl_linalg_balance_accum() */