lu.c 2.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128
  1. #include <math.h>
  2. #include <grass/gis.h>
  3. #include <grass/gmath.h>
  4. #define TINY 1.0e-20;
  5. /*!
  6. * \brief LU decomposition
  7. *
  8. * \param a double **
  9. * \param n int
  10. * \param indx int *
  11. * \param d double *
  12. *
  13. * \return 0 on singular matrix, 1 on success
  14. */
  15. int G_ludcmp(double **a, int n, int *indx, double *d)
  16. {
  17. int i, imax = 0, j, k;
  18. double big, dum, sum, temp;
  19. double *vv;
  20. int is_singular = FALSE;
  21. vv = G_alloc_vector(n);
  22. *d = 1.0;
  23. /* this pragma works, but doesn't really help speed things up */
  24. /* #pragma omp parallel for private(i, j, big, temp) shared(n, a, vv, is_singular) */
  25. for (i = 0; i < n; i++) {
  26. big = 0.0;
  27. for (j = 0; j < n; j++)
  28. if ((temp = fabs(a[i][j])) > big)
  29. big = temp;
  30. if (big == 0.0) {
  31. is_singular = TRUE;
  32. break;
  33. }
  34. vv[i] = 1.0 / big;
  35. }
  36. if (is_singular) {
  37. *d = 0.0;
  38. return 0; /* Singular matrix */
  39. }
  40. for (j = 0; j < n; j++) {
  41. for (i = 0; i < j; i++) {
  42. sum = a[i][j];
  43. for (k = 0; k < i; k++)
  44. sum -= a[i][k] * a[k][j];
  45. a[i][j] = sum;
  46. }
  47. big = 0.0;
  48. /* not very efficient, but this pragma helps speed things up a bit */
  49. #pragma omp parallel for private(i, k, sum, dum) shared(j, n, a, vv, big, imax)
  50. for (i = j; i < n; i++) {
  51. sum = a[i][j];
  52. for (k = 0; k < j; k++)
  53. sum -= a[i][k] * a[k][j];
  54. a[i][j] = sum;
  55. if ((dum = vv[i] * fabs(sum)) >= big) {
  56. big = dum;
  57. imax = i;
  58. }
  59. }
  60. if (j != imax) {
  61. for (k = 0; k < n; k++) {
  62. dum = a[imax][k];
  63. a[imax][k] = a[j][k];
  64. a[j][k] = dum;
  65. }
  66. *d = -(*d);
  67. vv[imax] = vv[j];
  68. }
  69. indx[j] = imax;
  70. if (a[j][j] == 0.0)
  71. a[j][j] = TINY;
  72. if (j != n) {
  73. dum = 1.0 / (a[j][j]);
  74. for (i = j + 1; i < n; i++)
  75. a[i][j] *= dum;
  76. }
  77. }
  78. G_free_vector(vv);
  79. return 1;
  80. }
  81. #undef TINY
  82. /*!
  83. * \brief LU backward substitution
  84. *
  85. * \param a double **
  86. * \param n int
  87. * \param indx int *
  88. * \param b double []
  89. *
  90. * \return void
  91. */
  92. void G_lubksb(double **a, int n, int *indx, double b[])
  93. {
  94. int i, ii, ip, j;
  95. double sum;
  96. ii = -1;
  97. for (i = 0; i < n; i++) {
  98. ip = indx[i];
  99. sum = b[ip];
  100. b[ip] = b[i];
  101. if (ii >= 0)
  102. for (j = ii; j < i; j++)
  103. sum -= a[i][j] * b[j];
  104. else if (sum)
  105. ii = i;
  106. b[i] = sum;
  107. }
  108. for (i = n - 1; i >= 0; i--) {
  109. sum = b[i];
  110. for (j = i + 1; j < n; j++)
  111. sum -= a[i][j] * b[j];
  112. b[i] = sum / a[i][i];
  113. }
  114. }