solvers_krylov.c 21 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807
  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. static G_math_spvector **create_diag_precond_matrix(double **A,
  25. G_math_spvector ** Asp,
  26. int rows, int prec);
  27. static int solver_pcg(double **A, G_math_spvector ** Asp, double *x,
  28. double *b, int rows, int maxit, double err, int prec, int has_band, int bandwidth);
  29. static int solver_cg(double **A, G_math_spvector ** Asp, double *x, double *b,
  30. int rows, int maxit, double err, int has_band, int bandwidth);
  31. static int solver_bicgstab(double **A, G_math_spvector ** Asp, double *x,
  32. double *b, int rows, int maxit, double err);
  33. /*!
  34. * \brief The iterative preconditioned conjugate gradients solver for symmetric positive definite matrices
  35. *
  36. * This iterative solver works with symmetric positive definite regular quadratic matrices.
  37. *
  38. * This solver solves the linear equation system:
  39. * A x = b
  40. *
  41. * The parameter <i>maxit</i> specifies the maximum number of iterations. If the maximum is reached, the
  42. * solver will abort the calculation and writes the current result into the vector x.
  43. * The parameter <i>err</i> defines the error break criteria for the solver.
  44. *
  45. * \param A (double **) -- the matrix
  46. * \param x (double *) -- the value vector
  47. * \param b (double *) -- the right hand side
  48. * \param rows (int)
  49. * \param maxit (int) -- the maximum number of iterations
  50. * \param err (double) -- defines the error break criteria
  51. * \param prec (int) -- the preconditioner which should be used 1,2 or 3
  52. * \return (int) -- 1 - success, 2 - not finished but success, 0 - matrix singular, -1 - could not solve the les
  53. *
  54. * */
  55. int G_math_solver_pcg(double **A, double *x, double *b, int rows, int maxit,
  56. double err, int prec)
  57. {
  58. return solver_pcg(A, NULL, x, b, rows, maxit, err, prec, 0, 0);
  59. }
  60. /*!
  61. * \brief The iterative preconditioned conjugate gradients solver for symmetric positive definite band matrices
  62. *
  63. * WARNING: The preconditioning of symmetric band matrices is not implemented yet
  64. *
  65. * This iterative solver works with symmetric positive definite band matrices.
  66. *
  67. * This solver solves the linear equation system:
  68. * A x = b
  69. *
  70. * The parameter <i>maxit</i> specifies the maximum number of iterations. If the maximum is reached, the
  71. * solver will abort the calculation and writes the current result into the vector x.
  72. * The parameter <i>err</i> defines the error break criteria for the solver.
  73. *
  74. * \param A (double **) -- the positive definite band matrix
  75. * \param x (double *) -- the value vector
  76. * \param b (double *) -- the right hand side
  77. * \param rows (int)
  78. * \param bandwidth (int) -- bandwidth of matrix A
  79. * \param maxit (int) -- the maximum number of iterations
  80. * \param err (double) -- defines the error break criteria
  81. * \param prec (int) -- the preconditioner which should be used 1,2 or 3
  82. * \return (int) -- 1 - success, 2 - not finished but success, 0 - matrix singular, -1 - could not solve the les
  83. *
  84. * */
  85. int G_math_solver_pcg_sband(double **A, double *x, double *b, int rows, int bandwidth, int maxit,
  86. double err, int prec)
  87. {
  88. G_fatal_error("Preconditioning of band matrics is not implemented yet");
  89. return solver_pcg(A, NULL, x, b, rows, maxit, err, prec, 1, bandwidth);
  90. }
  91. /*!
  92. * \brief The iterative preconditioned conjugate gradients solver for sparse symmetric positive definite matrices
  93. *
  94. * This iterative solver works with symmetric positive definite sparse matrices.
  95. *
  96. * This solver solves the linear equation system:
  97. * A x = b
  98. *
  99. * The parameter <i>maxit</i> specifies the maximum number of iterations. If the maximum is reached, the
  100. * solver will abort the calculation and writes the current result into the vector x.
  101. * The parameter <i>err</i> defines the error break criteria for the solver.
  102. *
  103. * \param Asp (G_math_spvector **) -- the sparse matrix
  104. * \param x (double *) -- the value vector
  105. * \param b (double *) -- the right hand side
  106. * \param rows (int)
  107. * \param maxit (int) -- the maximum number of iterations
  108. * \param err (double) -- defines the error break criteria
  109. * \param prec (int) -- the preconditioner which should be used 1,2 or 3
  110. * \return (int) -- 1 - success, 2 - not finished but success, 0 - matrix singular, -1 - could not solve the les
  111. *
  112. * */
  113. int G_math_solver_sparse_pcg(G_math_spvector ** Asp, double *x, double *b,
  114. int rows, int maxit, double err, int prec)
  115. {
  116. return solver_pcg(NULL, Asp, x, b, rows, maxit, err, prec, 0, 0);
  117. }
  118. int solver_pcg(double **A, G_math_spvector ** Asp, double *x, double *b,
  119. int rows, int maxit, double err, int prec, int has_band, int bandwidth)
  120. {
  121. double *r, *z;
  122. double *p;
  123. double *v;
  124. double s = 0.0;
  125. double a0 = 0, a1 = 0, mygamma, tmp = 0;
  126. int m, i;
  127. int finished = 2;
  128. int error_break;
  129. G_math_spvector **M;
  130. r = G_alloc_vector(rows);
  131. p = G_alloc_vector(rows);
  132. v = G_alloc_vector(rows);
  133. z = G_alloc_vector(rows);
  134. error_break = 0;
  135. /*compute the preconditioning matrix, this is a sparse matrix */
  136. M = create_diag_precond_matrix(A, Asp, rows, prec);
  137. /*
  138. * residual calculation
  139. */
  140. #pragma omp parallel
  141. {
  142. if (Asp)
  143. G_math_Ax_sparse(Asp, x, v, rows);
  144. else if(has_band)
  145. G_math_Ax_sband(A, x, v, rows, bandwidth);
  146. else
  147. G_math_d_Ax(A, x, v, rows, rows);
  148. G_math_d_ax_by(b, v, r, 1.0, -1.0, rows);
  149. /*performe the preconditioning */
  150. G_math_Ax_sparse(M, r, p, rows);
  151. /* scalar product */
  152. #pragma omp for schedule (static) private(i) reduction(+:s)
  153. for (i = 0; i < rows; i++) {
  154. s += p[i] * r[i];
  155. }
  156. }
  157. a0 = s;
  158. s = 0.0;
  159. /* ******************* */
  160. /* start the iteration */
  161. /* ******************* */
  162. for (m = 0; m < maxit; m++) {
  163. #pragma omp parallel default(shared)
  164. {
  165. if (Asp)
  166. G_math_Ax_sparse(Asp, p, v, rows);
  167. else if(has_band)
  168. G_math_Ax_sband(A, p, v, rows, bandwidth);
  169. else
  170. G_math_d_Ax(A, p, v, rows, rows);
  171. /* scalar product */
  172. #pragma omp for schedule (static) private(i) reduction(+:s)
  173. for (i = 0; i < rows; i++) {
  174. s += v[i] * p[i];
  175. }
  176. /* barrier */
  177. #pragma omp single
  178. {
  179. tmp = s;
  180. mygamma = a0 / tmp;
  181. s = 0.0;
  182. }
  183. G_math_d_ax_by(p, x, x, mygamma, 1.0, rows);
  184. if (m % 50 == 1) {
  185. if (Asp)
  186. G_math_Ax_sparse(Asp, x, v, rows);
  187. else if(has_band)
  188. G_math_Ax_sband(A, x, v, rows, bandwidth);
  189. else
  190. G_math_d_Ax(A, x, v, rows, rows);
  191. G_math_d_ax_by(b, v, r, 1.0, -1.0, rows);
  192. }
  193. else {
  194. G_math_d_ax_by(r, v, r, 1.0, -1.0 * mygamma, rows);
  195. }
  196. /*performe the preconditioning */
  197. G_math_Ax_sparse(M, r, z, rows);
  198. /* scalar product */
  199. #pragma omp for schedule (static) private(i) reduction(+:s)
  200. for (i = 0; i < rows; i++) {
  201. s += z[i] * r[i];
  202. }
  203. /* barrier */
  204. #pragma omp single
  205. {
  206. a1 = s;
  207. tmp = a1 / a0;
  208. a0 = a1;
  209. s = 0.0;
  210. if (a1 < 0 || a1 == 0 || a1 > 0) {
  211. ;
  212. }
  213. else {
  214. G_warning(_
  215. ("Unable to solve the linear equation system"));
  216. error_break = 1;
  217. }
  218. }
  219. G_math_d_ax_by(p, z, p, tmp, 1.0, rows);
  220. }
  221. if (Asp != NULL)
  222. G_message(_("Sparse PCG -- iteration %i error %g\n"), m, a0);
  223. else
  224. G_message(_("PCG -- iteration %i error %g\n"), m, a0);
  225. if (error_break == 1) {
  226. finished = -1;
  227. break;
  228. }
  229. if (a0 < err) {
  230. finished = 1;
  231. break;
  232. }
  233. }
  234. G_free(r);
  235. G_free(p);
  236. G_free(v);
  237. G_free(z);
  238. G_math_free_spmatrix(M, rows);
  239. return finished;
  240. }
  241. /*!
  242. * \brief The iterative conjugate gradients solver for symmetric positive definite matrices
  243. *
  244. * This iterative solver works with symmetric positive definite regular quadratic matrices.
  245. *
  246. * This solver solves the linear equation system:
  247. * A x = b
  248. *
  249. * The parameter <i>maxit</i> specifies the maximum number of iterations. If the maximum is reached, the
  250. * solver will abort the calculation and writes the current result into the vector x.
  251. * The parameter <i>err</i> defines the error break criteria for the solver.
  252. *
  253. * \param A (double **) -- the matrix
  254. * \param x (double *) -- the value vector
  255. * \param b (double *) -- the right hand side
  256. * \param rows (int)
  257. * \param maxit (int) -- the maximum number of iterations
  258. * \param err (double) -- defines the error break criteria
  259. * \return (int) -- 1 - success, 2 - not finished but success, 0 - matrix singular, -1 - could not solve the les
  260. *
  261. * */
  262. int G_math_solver_cg(double **A, double *x, double *b, int rows, int maxit,
  263. double err)
  264. {
  265. return solver_cg(A, NULL, x, b, rows, maxit, err, 0, 0);
  266. }
  267. /*!
  268. * \brief The iterative conjugate gradients solver for symmetric positive definite band matrices
  269. *
  270. * This iterative solver works with symmetric positive definite band matrices.
  271. *
  272. * This solver solves the linear equation system:
  273. * A x = b
  274. *
  275. * The parameter <i>maxit</i> specifies the maximum number of iterations. If the maximum is reached, the
  276. * solver will abort the calculation and writes the current result into the vector x.
  277. * The parameter <i>err</i> defines the error break criteria for the solver.
  278. *
  279. * \param A (double **) -- the symmetric positive definite band matrix
  280. * \param x (double *) -- the value vector
  281. * \param b (double *) -- the right hand side
  282. * \param rows (int)
  283. * \param bandwidth (int) -- the bandwidth of matrix A
  284. * \param maxit (int) -- the maximum number of iterations
  285. * \param err (double) -- defines the error break criteria
  286. * \return (int) -- 1 - success, 2 - not finished but success, 0 - matrix singular, -1 - could not solve the les
  287. *
  288. * */
  289. int G_math_solver_cg_sband(double **A, double *x, double *b, int rows, int bandwidth, int maxit, double err)
  290. {
  291. return solver_cg(A, NULL, x, b, rows, maxit, err, 1, bandwidth);
  292. }
  293. /*!
  294. * \brief The iterative conjugate gradients solver for sparse symmetric positive definite matrices
  295. *
  296. * This iterative solver works with symmetric positive definite sparse matrices.
  297. *
  298. * This solver solves the linear equation system:
  299. * A x = b
  300. *
  301. * The parameter <i>maxit</i> specifies the maximum number of iterations. If the maximum is reached, the
  302. * solver will abort the calculation and writes the current result into the vector x.
  303. * The parameter <i>err</i> defines the error break criteria for the solver.
  304. *
  305. * \param Asp (G_math_spvector **) -- the sparse matrix
  306. * \param x (double *) -- the value vector
  307. * \param b (double *) -- the right hand side
  308. * \param rows (int)
  309. * \param maxit (int) -- the maximum number of iterations
  310. * \param err (double) -- defines the error break criteria
  311. * \return (int) -- 1 - success, 2 - not finished but success, 0 - matrix singular, -1 - could not solve the les
  312. *
  313. * */
  314. int G_math_solver_sparse_cg(G_math_spvector ** Asp, double *x, double *b,
  315. int rows, int maxit, double err)
  316. {
  317. return solver_cg(NULL, Asp, x, b, rows, maxit, err, 0, 0);
  318. }
  319. int solver_cg(double **A, G_math_spvector ** Asp, double *x, double *b,
  320. int rows, int maxit, double err, int has_band, int bandwidth)
  321. {
  322. double *r;
  323. double *p;
  324. double *v;
  325. double s = 0.0;
  326. double a0 = 0, a1 = 0, mygamma, tmp = 0;
  327. int m, i;
  328. int finished = 2;
  329. int error_break;
  330. r = G_alloc_vector(rows);
  331. p = G_alloc_vector(rows);
  332. v = G_alloc_vector(rows);
  333. error_break = 0;
  334. /*
  335. * residual calculation
  336. */
  337. #pragma omp parallel
  338. {
  339. if (Asp)
  340. G_math_Ax_sparse(Asp, x, v, rows);
  341. else if(has_band)
  342. G_math_Ax_sband(A, x, v, rows, bandwidth);
  343. else
  344. G_math_d_Ax(A, x, v, rows, rows);
  345. G_math_d_ax_by(b, v, r, 1.0, -1.0, rows);
  346. G_math_d_copy(r, p, rows);
  347. /* scalar product */
  348. #pragma omp for schedule (static) private(i) reduction(+:s)
  349. for (i = 0; i < rows; i++) {
  350. s += r[i] * r[i];
  351. }
  352. }
  353. a0 = s;
  354. s = 0.0;
  355. /* ******************* */
  356. /* start the iteration */
  357. /* ******************* */
  358. for (m = 0; m < maxit; m++) {
  359. #pragma omp parallel default(shared)
  360. {
  361. if (Asp)
  362. G_math_Ax_sparse(Asp, p, v, rows);
  363. else if(has_band)
  364. G_math_Ax_sband(A, p, v, rows, bandwidth);
  365. else
  366. G_math_d_Ax(A, p, v, rows, rows);
  367. /* scalar product */
  368. #pragma omp for schedule (static) private(i) reduction(+:s)
  369. for (i = 0; i < rows; i++) {
  370. s += v[i] * p[i];
  371. }
  372. /* barrier */
  373. #pragma omp single
  374. {
  375. tmp = s;
  376. mygamma = a0 / tmp;
  377. s = 0.0;
  378. }
  379. G_math_d_ax_by(p, x, x, mygamma, 1.0, rows);
  380. if (m % 50 == 1) {
  381. if (Asp)
  382. G_math_Ax_sparse(Asp, x, v, rows);
  383. else if(has_band)
  384. G_math_Ax_sband(A, x, v, rows, bandwidth);
  385. else
  386. G_math_d_Ax(A, x, v, rows, rows);
  387. G_math_d_ax_by(b, v, r, 1.0, -1.0, rows);
  388. }
  389. else {
  390. G_math_d_ax_by(r, v, r, 1.0, -1.0 * mygamma, rows);
  391. }
  392. /* scalar product */
  393. #pragma omp for schedule (static) private(i) reduction(+:s)
  394. for (i = 0; i < rows; i++) {
  395. s += r[i] * r[i];
  396. }
  397. /* barrier */
  398. #pragma omp single
  399. {
  400. a1 = s;
  401. tmp = a1 / a0;
  402. a0 = a1;
  403. s = 0.0;
  404. if (a1 < 0 || a1 == 0 || a1 > 0) {
  405. ;
  406. }
  407. else {
  408. G_warning(_
  409. ("Unable to solve the linear equation system"));
  410. error_break = 1;
  411. }
  412. }
  413. G_math_d_ax_by(p, r, p, tmp, 1.0, rows);
  414. }
  415. if (Asp != NULL)
  416. G_message(_("Sparse CG -- iteration %i error %g\n"), m, a0);
  417. else
  418. G_message(_("CG -- iteration %i error %g\n"), m, a0);
  419. if (error_break == 1) {
  420. finished = -1;
  421. break;
  422. }
  423. if (a0 < err) {
  424. finished = 1;
  425. break;
  426. }
  427. }
  428. G_free(r);
  429. G_free(p);
  430. G_free(v);
  431. return finished;
  432. }
  433. /*!
  434. * \brief The iterative biconjugate gradients solver with stabilization for unsymmetric non-definite matrices
  435. *
  436. * This iterative solver works with regular quadratic matrices.
  437. *
  438. * This solver solves the linear equation system:
  439. * A x = b
  440. *
  441. * The parameter <i>maxit</i> specifies the maximum number of iterations. If the maximum is reached, the
  442. * solver will abort the calculation and writes the current result into the vector x.
  443. * The parameter <i>err</i> defines the error break criteria for the solver.
  444. *
  445. * \param A (double **) -- the matrix
  446. * \param x (double *) -- the value vector
  447. * \param b (double *) -- the right hand side
  448. * \param rows (int)
  449. * \param maxit (int) -- the maximum number of iterations
  450. * \param err (double) -- defines the error break criteria
  451. * \return (int) -- 1 - success, 2 - not finished but success, 0 - matrix singular, -1 - could not solve the les
  452. *
  453. * */
  454. int G_math_solver_bicgstab(double **A, double *x, double *b, int rows,
  455. int maxit, double err)
  456. {
  457. return solver_bicgstab(A, NULL, x, b, rows, maxit, err);
  458. }
  459. /*!
  460. * \brief The iterative biconjugate gradients solver with stabilization for unsymmetric non-definite matrices
  461. *
  462. * This iterative solver works with sparse matrices.
  463. *
  464. * This solver solves the linear equation system:
  465. * A x = b
  466. *
  467. * The parameter <i>maxit</i> specifies the maximum number of iterations. If the maximum is reached, the
  468. * solver will abort the calculation and writes the current result into the vector x.
  469. * The parameter <i>err</i> defines the error break criteria for the solver.
  470. *
  471. * \param Asp (G_math_spvector **) -- the sparse matrix
  472. * \param x (double *) -- the value vector
  473. * \param b (double *) -- the right hand side
  474. * \param rows (int)
  475. * \param maxit (int) -- the maximum number of iterations
  476. * \param err (double) -- defines the error break criteria
  477. * \return (int) -- 1 - success, 2 - not finished but success, 0 - matrix singular, -1 - could not solve the les
  478. *
  479. * */
  480. int G_math_solver_sparse_bicgstab(G_math_spvector ** Asp, double *x,
  481. double *b, int rows, int maxit, double err)
  482. {
  483. return solver_bicgstab(NULL, Asp, x, b, rows, maxit, err);
  484. }
  485. int solver_bicgstab(double **A, G_math_spvector ** Asp, double *x, double *b,
  486. int rows, int maxit, double err)
  487. {
  488. double *r;
  489. double *r0;
  490. double *p;
  491. double *v;
  492. double *s;
  493. double *t;
  494. double s1 = 0.0, s2 = 0.0, s3 = 0.0;
  495. double alpha = 0, beta = 0, omega, rr0 = 0, error;
  496. int m, i;
  497. int finished = 2;
  498. int error_break;
  499. r = G_alloc_vector(rows);
  500. r0 = G_alloc_vector(rows);
  501. p = G_alloc_vector(rows);
  502. v = G_alloc_vector(rows);
  503. s = G_alloc_vector(rows);
  504. t = G_alloc_vector(rows);
  505. error_break = 0;
  506. #pragma omp parallel
  507. {
  508. if (Asp)
  509. G_math_Ax_sparse(Asp, x, v, rows);
  510. else
  511. G_math_d_Ax(A, x, v, rows, rows);
  512. G_math_d_ax_by(b, v, r, 1.0, -1.0, rows);
  513. G_math_d_copy(r, r0, rows);
  514. G_math_d_copy(r, p, rows);
  515. }
  516. s1 = s2 = s3 = 0.0;
  517. /* ******************* */
  518. /* start the iteration */
  519. /* ******************* */
  520. for (m = 0; m < maxit; m++) {
  521. #pragma omp parallel default(shared)
  522. {
  523. if (Asp)
  524. G_math_Ax_sparse(Asp, p, v, rows);
  525. else
  526. G_math_d_Ax(A, p, v, rows, rows);
  527. /* scalar product */
  528. #pragma omp for schedule (static) private(i) reduction(+:s1, s2, s3)
  529. for (i = 0; i < rows; i++) {
  530. s1 += r[i] * r[i];
  531. s2 += r[i] * r0[i];
  532. s3 += v[i] * r0[i];
  533. }
  534. #pragma omp single
  535. {
  536. error = s1;
  537. if (error < 0 || error == 0 || error > 0) {
  538. ;
  539. }
  540. else {
  541. G_warning(_
  542. ("Unable to solve the linear equation system"));
  543. error_break = 1;
  544. }
  545. rr0 = s2;
  546. alpha = rr0 / s3;
  547. s1 = s2 = s3 = 0.0;
  548. }
  549. G_math_d_ax_by(r, v, s, 1.0, -1.0 * alpha, rows);
  550. if (Asp)
  551. G_math_Ax_sparse(Asp, s, t, rows);
  552. else
  553. G_math_d_Ax(A, s, t, rows, rows);
  554. /* scalar product */
  555. #pragma omp for schedule (static) private(i) reduction(+:s1, s2)
  556. for (i = 0; i < rows; i++) {
  557. s1 += t[i] * s[i];
  558. s2 += t[i] * t[i];
  559. }
  560. #pragma omp single
  561. {
  562. omega = s1 / s2;
  563. s1 = s2 = 0.0;
  564. }
  565. G_math_d_ax_by(p, s, r, alpha, omega, rows);
  566. G_math_d_ax_by(x, r, x, 1.0, 1.0, rows);
  567. G_math_d_ax_by(s, t, r, 1.0, -1.0 * omega, rows);
  568. #pragma omp for schedule (static) private(i) reduction(+:s1)
  569. for (i = 0; i < rows; i++) {
  570. s1 += r[i] * r0[i];
  571. }
  572. #pragma omp single
  573. {
  574. beta = alpha / omega * s1 / rr0;
  575. s1 = s2 = s3 = 0.0;
  576. }
  577. G_math_d_ax_by(p, v, p, 1.0, -1.0 * omega, rows);
  578. G_math_d_ax_by(p, r, p, beta, 1.0, rows);
  579. }
  580. if (Asp != NULL)
  581. G_message(_("Sparse BiCGStab -- iteration %i error %g\n"), m,
  582. error);
  583. else
  584. G_message(_("BiCGStab -- iteration %i error %g\n"), m, error);
  585. if (error_break == 1) {
  586. finished = -1;
  587. break;
  588. }
  589. if (error < err) {
  590. finished = 1;
  591. break;
  592. }
  593. }
  594. G_free(r);
  595. G_free(r0);
  596. G_free(p);
  597. G_free(v);
  598. G_free(s);
  599. G_free(t);
  600. return finished;
  601. }
  602. /*!
  603. * \brief Compute a diagonal preconditioning matrix for krylov space solver
  604. *
  605. * \param A (double **) -- the matrix for which the precondition should be computed (if the sparse matrix is used, set it to NULL)
  606. * \param Asp (G_math_spvector **) -- the matrix for which the precondition should be computed
  607. * \param rows (int)
  608. * \param prec (int) -- which preconditioner should be used 1, 2 or 3
  609. *
  610. * */
  611. G_math_spvector **create_diag_precond_matrix(double **A,
  612. G_math_spvector ** Asp, int rows,
  613. int prec)
  614. {
  615. G_math_spvector **Msp;
  616. int i, j, cols = rows;
  617. double sum;
  618. Msp = G_math_alloc_spmatrix(rows);
  619. if (A != NULL) {
  620. #pragma omp parallel for schedule (static) private(i, j, sum) shared(A, Msp, rows, cols, prec)
  621. for (i = 0; i < rows; i++) {
  622. G_math_spvector *spvect = G_math_alloc_spvector(1);
  623. switch (prec) {
  624. case G_MATH_ROWSCALE_EUKLIDNORM_PRECONDITION:
  625. sum = 0;
  626. for (j = 0; j < cols; j++)
  627. sum += A[i][j] * A[i][j];
  628. spvect->values[0] = 1.0 / sqrt(sum);
  629. break;
  630. case G_MATH_ROWSCALE_ABSSUMNORM_PRECONDITION:
  631. sum = 0;
  632. for (j = 0; j < cols; j++)
  633. sum += fabs(A[i][j]);
  634. spvect->values[0] = 1.0 / (sum);
  635. break;
  636. case G_MATH_DIAGONAL_PRECONDITION:
  637. default:
  638. spvect->values[0] = 1.0 / A[i][i];
  639. break;
  640. }
  641. spvect->index[0] = i;
  642. spvect->cols = 1;;
  643. G_math_add_spvector(Msp, spvect, i);
  644. }
  645. }
  646. else {
  647. #pragma omp parallel for schedule (static) private(i, j, sum) shared(Asp, Msp, rows, cols, prec)
  648. for (i = 0; i < rows; i++) {
  649. G_math_spvector *spvect = G_math_alloc_spvector(1);
  650. switch (prec) {
  651. case G_MATH_ROWSCALE_EUKLIDNORM_PRECONDITION:
  652. sum = 0;
  653. for (j = 0; j < Asp[i]->cols; j++)
  654. sum += Asp[i]->values[j] * Asp[i]->values[j];
  655. spvect->values[0] = 1.0 / sqrt(sum);
  656. break;
  657. case G_MATH_ROWSCALE_ABSSUMNORM_PRECONDITION:
  658. sum = 0;
  659. for (j = 0; j < Asp[i]->cols; j++)
  660. sum += fabs(Asp[i]->values[j]);
  661. spvect->values[0] = 1.0 / (sum);
  662. break;
  663. case G_MATH_DIAGONAL_PRECONDITION:
  664. default:
  665. for (j = 0; j < Asp[i]->cols; j++)
  666. if (i == Asp[i]->index[j])
  667. spvect->values[0] = 1.0 / Asp[i]->values[j];
  668. break;
  669. }
  670. spvect->index[0] = i;
  671. spvect->cols = 1;;
  672. G_math_add_spvector(Msp, spvect, i);
  673. }
  674. }
  675. return Msp;
  676. }