solvers_direct.c 7.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329
  1. /*****************************************************************************
  2. *
  3. * MODULE: Grass numerical math interface
  4. * AUTHOR(S): Soeren Gebbert, Berlin (GER) Dec 2006
  5. * soerengebbert <at> googlemail <dot> com
  6. *
  7. * PURPOSE: linear equation system solvers
  8. * part of the gmath library
  9. *
  10. * COPYRIGHT: (C) 2010 by the GRASS Development Team
  11. *
  12. * This program is free software under the GNU General Public
  13. * License (>=v2). Read the file COPYING that comes with GRASS
  14. * for details.
  15. *
  16. *****************************************************************************/
  17. #include <math.h>
  18. #include <unistd.h>
  19. #include <stdio.h>
  20. #include <string.h>
  21. #include <grass/gis.h>
  22. #include <grass/gmath.h>
  23. #include <grass/glocale.h>
  24. #define TINY 1.0e-20
  25. #define COMP_PIVOT 100
  26. /*!
  27. * \brief The gauss elimination solver for quardatic matrices
  28. *
  29. * This solver does not support sparse matrices
  30. * The matrix A will be overwritten.
  31. * The result is written to the vector x
  32. *
  33. * \param A double **
  34. * \param x double *
  35. * \param b double *
  36. * \param rows int
  37. * \return int -- 1 success
  38. * */
  39. int G_math_solver_gauss(double **A, double *x, double *b, int rows)
  40. {
  41. G_message(_("Starting direct gauss elimination solver"));
  42. G_math_gauss_elimination(A, b, rows);
  43. G_math_backward_substitution(A, x, b, rows);
  44. return 1;
  45. }
  46. /*!
  47. * \brief The LU solver for quardatic matrices
  48. *
  49. * This solver does not support sparse matrices
  50. * The matrix A will be overwritten.
  51. * The result is written to the vector x in the G_math_les structure
  52. *
  53. *
  54. * \param A double **
  55. * \param x double *
  56. * \param b double *
  57. * \param rows int
  58. * \return int -- 1 success
  59. * */
  60. int G_math_solver_lu(double **A, double *x, double *b, int rows)
  61. {
  62. int i;
  63. double *c, *tmpv;
  64. G_message(_("Starting direct lu decomposition solver"));
  65. tmpv = G_alloc_vector(rows);
  66. c = G_alloc_vector(rows);
  67. G_math_lu_decomposition(A, b, rows);
  68. #pragma omp parallel
  69. {
  70. #pragma omp for schedule (static) private(i)
  71. for (i = 0; i < rows; i++) {
  72. tmpv[i] = A[i][i];
  73. A[i][i] = 1;
  74. }
  75. #pragma omp single
  76. {
  77. G_math_forward_substitution(A, b, b, rows);
  78. }
  79. #pragma omp for schedule (static) private(i)
  80. for (i = 0; i < rows; i++) {
  81. A[i][i] = tmpv[i];
  82. }
  83. #pragma omp single
  84. {
  85. G_math_backward_substitution(A, x, b, rows);
  86. }
  87. }
  88. G_free(c);
  89. G_free(tmpv);
  90. return 1;
  91. }
  92. /*!
  93. * \brief The choleksy decomposition solver for quardatic, symmetric
  94. * positiv definite matrices
  95. *
  96. * This solver does not support sparse matrices
  97. * The matrix A will be overwritten.
  98. * The result is written to the vector x
  99. *
  100. * \param A double **
  101. * \param x double *
  102. * \param b double *
  103. * \param bandwith int -- the bandwith of the band matrix, if unsure set to rows
  104. * \param rows int
  105. * \return int -- 1 success
  106. * */
  107. int G_math_solver_cholesky(double **A, double *x, double *b, int bandwith,
  108. int rows)
  109. {
  110. G_message(_("Starting cholesky decomposition solver"));
  111. if (G_math_cholesky_decomposition(A, rows, bandwith) != 1) {
  112. G_warning(_("Unable to solve the linear equation system"));
  113. return -2;
  114. }
  115. G_math_forward_substitution(A, b, b, rows);
  116. G_math_backward_substitution(A, x, b, rows);
  117. return 1;
  118. }
  119. /*!
  120. * \brief Gauss elimination
  121. *
  122. * To run this solver efficiently,
  123. * no pivoting is supported.
  124. * The matrix will be overwritten with the decomposite form
  125. * \param A double **
  126. * \param b double *
  127. * \param rows int
  128. * \return void
  129. *
  130. * */
  131. void G_math_gauss_elimination(double **A, double *b, int rows)
  132. {
  133. int i, j, k;
  134. double tmpval = 0.0;
  135. for (k = 0; k < rows - 1; k++) {
  136. #pragma omp parallel for schedule (static) private(i, j, tmpval) shared(k, A, b, rows)
  137. for (i = k + 1; i < rows; i++) {
  138. tmpval = A[i][k] / A[k][k];
  139. b[i] = b[i] - tmpval * b[k];
  140. for (j = k + 1; j < rows; j++) {
  141. A[i][j] = A[i][j] - tmpval * A[k][j];
  142. }
  143. }
  144. }
  145. return;
  146. }
  147. /*!
  148. * \brief lu decomposition
  149. *
  150. * To run this solver efficiently,
  151. * no pivoting is supported.
  152. * The matrix will be overwritten with the decomposite form
  153. *
  154. * \param A double **
  155. * \param b double * -- this vector is needed if its part of the linear equation system, otherwise set it to NULL
  156. * \param rows int
  157. * \return void
  158. *
  159. * */
  160. void G_math_lu_decomposition(double **A, double *b, int rows)
  161. {
  162. int i, j, k;
  163. for (k = 0; k < rows - 1; k++) {
  164. #pragma omp parallel for schedule (static) private(i, j) shared(k, A, rows)
  165. for (i = k + 1; i < rows; i++) {
  166. A[i][k] = A[i][k] / A[k][k];
  167. for (j = k + 1; j < rows; j++) {
  168. A[i][j] = A[i][j] - A[i][k] * A[k][j];
  169. }
  170. }
  171. }
  172. return;
  173. }
  174. /*!
  175. * \brief cholesky decomposition for symmetric, positiv definite matrices
  176. * with bandwith optimization
  177. *
  178. * The provided matrix will be overwritten with the lower and
  179. * upper triangle matrix A = LL^T
  180. *
  181. * \param A double **
  182. * \param rows int
  183. * \param bandwith int -- the bandwith of the matrix (0 > bandwith <= cols)
  184. * \return void
  185. *
  186. * */
  187. int G_math_cholesky_decomposition(double **A, int rows, int bandwith)
  188. {
  189. int i = 0, j = 0, k = 0;
  190. double sum_1 = 0.0;
  191. double sum_2 = 0.0;
  192. int colsize;
  193. if (bandwith <= 0)
  194. bandwith = rows;
  195. colsize = bandwith;
  196. for (k = 0; k < rows; k++) {
  197. #pragma omp parallel for schedule (static) private(i, j, sum_2) shared(A, k) reduction(+:sum_1)
  198. for (j = 0; j < k; j++) {
  199. sum_1 += A[k][j] * A[k][j];
  200. }
  201. if (0 > (A[k][k] - sum_1)) {
  202. G_warning("Matrix is not positive definite. break.");
  203. return -1;
  204. }
  205. A[k][k] = sqrt(A[k][k] - sum_1);
  206. sum_1 = 0.0;
  207. if ((k + bandwith) > rows) {
  208. colsize = rows;
  209. }
  210. else {
  211. colsize = k + bandwith;
  212. }
  213. #pragma omp parallel for schedule (static) private(i, j, sum_2) shared(A, k, sum_1, colsize)
  214. for (i = k + 1; i < colsize; i++) {
  215. sum_2 = 0.0;
  216. for (j = 0; j < k; j++) {
  217. sum_2 += A[i][j] * A[k][j];
  218. }
  219. A[i][k] = (A[i][k] - sum_2) / A[k][k];
  220. }
  221. }
  222. /* we need to copy the lower triangle matrix to the upper triangle */
  223. #pragma omp parallel for schedule (static) private(i, k) shared(A, rows)
  224. for (k = 0; k < rows; k++) {
  225. for (i = k + 1; i < rows; i++) {
  226. A[k][i] = A[i][k];
  227. }
  228. }
  229. return 1;
  230. }
  231. /*!
  232. * \brief backward substitution
  233. *
  234. * \param A double **
  235. * \param x double *
  236. * \param b double *
  237. * \param rows int
  238. * \return void
  239. *
  240. * */
  241. void G_math_backward_substitution(double **A, double *x, double *b, int rows)
  242. {
  243. int i, j;
  244. for (i = rows - 1; i >= 0; i--) {
  245. for (j = i + 1; j < rows; j++) {
  246. b[i] = b[i] - A[i][j] * x[j];
  247. }
  248. x[i] = (b[i]) / A[i][i];
  249. }
  250. return;
  251. }
  252. /*!
  253. * \brief forward substitution
  254. *
  255. * \param A double **
  256. * \param x double *
  257. * \param b double *
  258. * \param rows int
  259. * \return void
  260. *
  261. * */
  262. void G_math_forward_substitution(double **A, double *x, double *b, int rows)
  263. {
  264. int i, j;
  265. double tmpval = 0.0;
  266. for (i = 0; i < rows; i++) {
  267. tmpval = 0;
  268. for (j = 0; j < i; j++) {
  269. tmpval += A[i][j] * x[j];
  270. }
  271. x[i] = (b[i] - tmpval) / A[i][i];
  272. }
  273. return;
  274. }