solvers_direct.c 9.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417
  1. /*****************************************************************************
  2. *
  3. * MODULE: Grass PDE Numerical Library
  4. * AUTHOR(S): Soeren Gebbert, Berlin (GER) Dec 2006
  5. * soerengebbert <at> gmx <dot> de
  6. *
  7. * PURPOSE: direkt linear equation system solvers
  8. * part of the gpde library
  9. *
  10. * COPYRIGHT: (C) 2007 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/glocale.h"
  23. #include "grass/gmath.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. * \int 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_solving(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. * \int 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_solving(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_solving(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. * \int rows int
  104. * \return int -- 1 success
  105. * */
  106. int G_math_solver_cholesky(double **A, double *x, double *b, int bandwith,
  107. int rows)
  108. {
  109. G_message(_("Starting cholesky decomposition solver"));
  110. if (G_math_cholesky_decomposition(A, rows, bandwith) != 1) {
  111. G_warning(_("Unable to solve the linear equation system"));
  112. return -2;
  113. }
  114. G_math_forward_solving(A, b, b, rows);
  115. G_math_backward_solving(A, x, b, rows);
  116. return 1;
  117. }
  118. /*!
  119. * \brief Gauss elimination
  120. *
  121. * To run this solver efficiently,
  122. * no pivoting is supported.
  123. * The matrix will be overwritten with the decomposite form
  124. * \param A double **
  125. * \param b double *
  126. * \param rows int
  127. * \return void
  128. *
  129. * */
  130. void G_math_gauss_elimination(double **A, double *b, int rows)
  131. {
  132. int i, j, k;
  133. double tmpval = 0.0;
  134. /*compute the pivot -- commented out, because its meaningless
  135. to compute it only nth times. */
  136. /*G_math_pivot_create(A, b, rows, 0); */
  137. for (k = 0; k < rows - 1; k++) {
  138. #pragma omp parallel for schedule (static) private(i, j, tmpval) shared(k, A, b, rows)
  139. for (i = k + 1; i < rows; i++) {
  140. tmpval = A[i][k] / A[k][k];
  141. b[i] = b[i] - tmpval * b[k];
  142. for (j = k + 1; j < rows; j++) {
  143. A[i][j] = A[i][j] - tmpval * A[k][j];
  144. }
  145. }
  146. }
  147. return;
  148. }
  149. /*!
  150. * \brief lu decomposition
  151. *
  152. * To run this solver efficiently,
  153. * no pivoting is supported.
  154. * The matrix will be overwritten with the decomposite form
  155. *
  156. * \param A double **
  157. * \param b double * -- this vector is needed if its part of the linear equation system, otherwise set it to NULL
  158. * \param rows int
  159. * \return void
  160. *
  161. * */
  162. void G_math_lu_decomposition(double **A, double *b, int rows)
  163. {
  164. int i, j, k;
  165. /*compute the pivot -- commented out, because its meaningless
  166. to compute it only nth times. */
  167. /*G_math_pivot_create(A, b, rows, 0); */
  168. for (k = 0; k < rows - 1; k++) {
  169. #pragma omp parallel for schedule (static) private(i, j) shared(k, A, rows)
  170. for (i = k + 1; i < rows; i++) {
  171. A[i][k] = A[i][k] / A[k][k];
  172. for (j = k + 1; j < rows; j++) {
  173. A[i][j] = A[i][j] - A[i][k] * A[k][j];
  174. }
  175. }
  176. }
  177. return;
  178. }
  179. /*!
  180. * \brief cholesky decomposition for symmetric, positiv definite matrices
  181. * with bandwith optimization
  182. *
  183. * The provided matrix will be overwritten with the lower and
  184. * upper triangle matrix A = LL^T
  185. *
  186. * \param A double **
  187. * \param rows int
  188. * \param bandwith int -- the bandwith of the matrix (0 > bandwith <= cols)
  189. * \return void
  190. *
  191. * */
  192. int G_math_cholesky_decomposition(double **A, int rows, int bandwith)
  193. {
  194. int i = 0, j = 0, k = 0;
  195. double sum_1 = 0.0;
  196. double sum_2 = 0.0;
  197. int colsize;
  198. if (bandwith <= 0)
  199. bandwith = rows;
  200. colsize = bandwith;
  201. for (k = 0; k < rows; k++) {
  202. #pragma omp parallel for schedule (static) private(i, j, sum_2) shared(A, k) reduction(+:sum_1)
  203. for (j = 0; j < k; j++) {
  204. sum_1 += A[k][j] * A[k][j];
  205. }
  206. if (0 > (A[k][k] - sum_1)) {
  207. G_warning("Matrix is not positive definite. break.");
  208. return -1;
  209. }
  210. A[k][k] = sqrt(A[k][k] - sum_1);
  211. sum_1 = 0.0;
  212. if ((k + bandwith) > rows) {
  213. colsize = rows;
  214. }
  215. else {
  216. colsize = k + bandwith;
  217. }
  218. #pragma omp parallel for schedule (static) private(i, j, sum_2) shared(A, k, sum_1, colsize)
  219. for (i = k + 1; i < colsize; i++) {
  220. sum_2 = 0.0;
  221. for (j = 0; j < k; j++) {
  222. sum_2 += A[i][j] * A[k][j];
  223. }
  224. A[i][k] = (A[i][k] - sum_2) / A[k][k];
  225. }
  226. }
  227. /*we need to copy the lower triangle matrix to the upper trianle */
  228. #pragma omp parallel for schedule (static) private(i, k) shared(A, rows)
  229. for (k = 0; k < rows; k++) {
  230. for (i = k + 1; i < rows; i++) {
  231. A[k][i] = A[i][k];
  232. }
  233. }
  234. return 1;
  235. }
  236. /*!
  237. * \brief backward solving
  238. *
  239. * \param A double **
  240. * \param x double *
  241. * \param b double *
  242. * \param rows int
  243. * \return void
  244. *
  245. * */
  246. void G_math_backward_solving(double **A, double *x, double *b, int rows)
  247. {
  248. int i, j;
  249. for (i = rows - 1; i >= 0; i--) {
  250. for (j = i + 1; j < rows; j++) {
  251. b[i] = b[i] - A[i][j] * x[j];
  252. }
  253. x[i] = (b[i]) / A[i][i];
  254. }
  255. return;
  256. }
  257. /*!
  258. * \brief forward solving
  259. *
  260. * \param A double **
  261. * \param x double *
  262. * \param b double *
  263. * \param rows int
  264. * \return void
  265. *
  266. * */
  267. void G_math_forward_solving(double **A, double *x, double *b, int rows)
  268. {
  269. int i, j;
  270. double tmpval = 0.0;
  271. for (i = 0; i < rows; i++) {
  272. tmpval = 0;
  273. for (j = 0; j < i; j++) {
  274. tmpval += A[i][j] * x[j];
  275. }
  276. x[i] = (b[i] - tmpval) / A[i][i];
  277. }
  278. return;
  279. }
  280. /*!
  281. * \brief Optimize the structure of the linear equation system with a common pivoting strategy
  282. *
  283. * Create a optimized linear equation system for
  284. * direct solvers: gauss and lu decomposition.
  285. *
  286. * The rows are permuted based on the pivot elements.
  287. *
  288. * This algorithm will modify the provided linear equation system
  289. * and should only be used with the gauss elimination and lu decomposition solver.
  290. *
  291. * \param A double ** - a quadratic matrix
  292. * \param b double * - the right hand vector, if not available set it to NULL
  293. * \param rows int
  294. * \param start int -- the row
  295. * \return int - the number of swapped rows
  296. *
  297. *
  298. * */
  299. int G_math_pivot_create(double **A, double *b, int rows, int start)
  300. {
  301. int num = 0; /*number of changed rows */
  302. int i, j, k;
  303. double max;
  304. int number = 0;
  305. double tmpval = 0.0, s = 0.0;
  306. double *link = NULL;
  307. link = G_alloc_vector(rows);
  308. G_debug(2, "G_math_pivot_create: swap rows if needed");
  309. for (i = start; i < rows; i++) {
  310. s = 0.0;
  311. for (k = i + 1; k < rows; k++) {
  312. s += fabs(A[i][k]);
  313. }
  314. max = fabs(A[i][i]) / s;
  315. number = i;
  316. for (j = i + 1; j < rows; j++) {
  317. s = 0.0;
  318. for (k = j; k < rows; k++) {
  319. s += fabs(A[j][k]);
  320. }
  321. /*search for the pivot element */
  322. if (max < fabs(A[j][i]) / s) {
  323. max = fabs(A[j][i] / s);
  324. number = j;
  325. }
  326. }
  327. if (max == 0) {
  328. max = TINY;
  329. G_warning("Matrix is singular");
  330. }
  331. /*if an pivot element was found, swap the les entries */
  332. if (number != i) {
  333. G_debug(4, "swap row %i with row %i", i, number);
  334. if (b != NULL) {
  335. tmpval = b[number];
  336. b[number] = b[i];
  337. b[i] = tmpval;
  338. }
  339. G_math_d_copy(A[number], link, rows);
  340. G_math_d_copy(A[i], A[number], rows);
  341. G_math_d_copy(link, A[i], rows);
  342. num++;
  343. }
  344. }
  345. G_free_vector(link);
  346. return num;
  347. }