test_solvers.c 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507
  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: Unit tests for les solving
  8. *
  9. * COPYRIGHT: (C) 2000 by the GRASS Development Team
  10. *
  11. * This program is free software under the GNU General Public
  12. * License (>=v2). Read the file COPYING that comes with GRASS
  13. * for details.
  14. *
  15. *****************************************************************************/
  16. #include <stdio.h>
  17. #include <stdlib.h>
  18. #include <string.h>
  19. #include <grass/glocale.h>
  20. #include <grass/gmath.h>
  21. #include "test_gmath_lib.h"
  22. #define EPSILON_DIRECT 1.0E-10
  23. #define EPSILON_ITER 1.0E-4
  24. /* prototypes */
  25. static int test_solvers(void);
  26. /* ************************************************************************* */
  27. /* Performe the solver unit tests ****************************************** */
  28. /* ************************************************************************* */
  29. int unit_test_solvers(void)
  30. {
  31. int sum = 0;
  32. G_message(_("\n++ Running solver unit tests ++"));
  33. sum += test_solvers();
  34. if (sum > 0)
  35. G_warning(_("\n-- Solver unit tests failure --"));
  36. else
  37. G_message(_("\n-- Solver unit tests finished successfully --"));
  38. return sum;
  39. }
  40. /* *************************************************************** */
  41. /* Test all implemented solvers for sparse and normal matrix *** */
  42. /* *************************************************************** */
  43. int test_solvers(void)
  44. {
  45. G_math_les *les;
  46. G_math_les *sples;
  47. int sum = 0;
  48. double val = 0.0;
  49. G_message("\t * testing jacobi solver with symmetric matrix\n");
  50. les = create_normal_symmetric_les(TEST_NUM_ROWS);
  51. sples = create_sparse_symmetric_les(TEST_NUM_ROWS);
  52. G_math_solver_jacobi(les->A, les->x, les->b, les->rows, 250, 1, 0.1e-10);
  53. G_math_d_asum_norm(les->x, &val, les->rows);
  54. if ((val - (double)les->rows) > EPSILON_ITER)
  55. {
  56. G_warning("Error in G_math_solver_jacobi abs %2.20f != %i", val,
  57. les->rows);
  58. sum++;
  59. }
  60. G_math_solver_sparse_jacobi(sples->Asp, sples->x, sples->b, les->rows, 250,
  61. 1, 0.1e-10);
  62. G_math_d_asum_norm(sples->x, &val, sples->rows);
  63. if ((val - (double)sples->rows) > EPSILON_ITER)
  64. {
  65. G_warning("Error in G_math_solver_sparse_jacobi abs %2.20f != %i", val,
  66. sples->rows);
  67. sum++;
  68. }
  69. G_math_free_les(les);
  70. G_math_free_les(sples);
  71. G_message("\t * testing jacobi solver with unsymmetric matrix\n");
  72. les = create_normal_unsymmetric_les(TEST_NUM_ROWS);
  73. sples = create_sparse_unsymmetric_les(TEST_NUM_ROWS);
  74. G_math_solver_jacobi(les->A, les->x, les->b, les->rows, 250, 1, 0.1e-10);
  75. G_math_d_asum_norm(les->x, &val, les->rows);
  76. if ((val - (double)les->rows) > EPSILON_ITER)
  77. {
  78. G_warning("Error in G_math_solver_jacobi abs %2.20f != %i", val,
  79. les->rows);
  80. sum++;
  81. }
  82. G_math_solver_sparse_jacobi(sples->Asp, sples->x, sples->b, les->rows, 250,
  83. 1, 0.1e-10);
  84. G_math_d_asum_norm(sples->x, &val, sples->rows);
  85. if ((val - (double)sples->rows) > EPSILON_ITER)
  86. {
  87. G_warning("Error in G_math_solver_sparse_jacobi abs %2.20f != %i", val,
  88. sples->rows);
  89. sum++;
  90. }
  91. G_math_free_les(les);
  92. G_math_free_les(sples);
  93. G_message("\t * testing gauss seidel solver with symmetric matrix\n");
  94. les = create_normal_symmetric_les(TEST_NUM_ROWS);
  95. sples = create_sparse_symmetric_les(TEST_NUM_ROWS);
  96. G_math_solver_gs(les->A, les->x, les->b, les->rows, 150, 1, 0.1e-9);
  97. G_math_d_asum_norm(les->x, &val, les->rows);
  98. if ((val - (double)les->rows) > EPSILON_ITER)
  99. {
  100. G_warning("Error in G_math_solver_gs abs %2.20f != %i", val, les->rows);
  101. sum++;
  102. }
  103. G_math_solver_sparse_gs(sples->Asp, sples->x, sples->b, les->rows, 150, 1,
  104. 0.1e-9);
  105. G_math_d_asum_norm(sples->x, &val, sples->rows);
  106. if ((val - (double)sples->rows) > EPSILON_ITER)
  107. {
  108. G_warning("Error in G_math_solver_sparse_gs abs %2.20f != %i", val,
  109. sples->rows);
  110. sum++;
  111. }
  112. G_math_free_les(les);
  113. G_math_free_les(sples);
  114. G_message("\t * testing gauss seidel solver with unsymmetric matrix\n");
  115. les = create_normal_unsymmetric_les(TEST_NUM_ROWS);
  116. sples = create_sparse_unsymmetric_les(TEST_NUM_ROWS);
  117. G_math_solver_gs(les->A, les->x, les->b, les->rows, 150, 1, 0.1e-9);
  118. G_math_d_asum_norm(les->x, &val, les->rows);
  119. if ((val - (double)les->rows) > EPSILON_ITER)
  120. {
  121. G_warning("Error in G_math_solver_gs abs %2.20f != %i", val, les->rows);
  122. sum++;
  123. }
  124. G_math_solver_sparse_gs(sples->Asp, sples->x, sples->b, les->rows, 150, 1,
  125. 0.1e-9);
  126. G_math_d_asum_norm(sples->x, &val, sples->rows);
  127. if ((val - (double)sples->rows) > EPSILON_ITER)
  128. {
  129. G_warning("Error in G_math_solver_sparse_gs abs %2.20f != %i", val,
  130. sples->rows);
  131. sum++;
  132. }
  133. G_math_free_les(les);
  134. G_math_free_les(sples);
  135. G_message("\t * testing pcg solver with symmetric bad conditioned matrix and preconditioner 3\n");
  136. les = create_normal_symmetric_pivot_les(TEST_NUM_ROWS);
  137. G_math_solver_pcg(les->A, les->x, les->b, les->rows, 250, 0.1e-9, 3);
  138. G_math_d_asum_norm(les->x, &val, les->rows);
  139. if ((val - (double)les->rows) > EPSILON_ITER)
  140. {
  141. G_warning("Error in G_math_solver_pcg abs %2.20f != %i", val, les->rows);
  142. sum++;
  143. }
  144. G_math_print_les(les);
  145. G_math_free_les(les);
  146. G_message("\t * testing pcg solver with symmetric matrix and preconditioner 1\n");
  147. les = create_normal_symmetric_les(TEST_NUM_ROWS);
  148. sples = create_sparse_symmetric_les(TEST_NUM_ROWS);
  149. G_math_solver_pcg(les->A, les->x, les->b, les->rows, 250, 0.1e-9, 1);
  150. G_math_d_asum_norm(les->x, &val, les->rows);
  151. if ((val - (double)les->rows) > EPSILON_ITER)
  152. {
  153. G_warning("Error in G_math_solver_pcg abs %2.20f != %i", val, les->rows);
  154. sum++;
  155. }
  156. G_math_print_les(les);
  157. G_math_solver_sparse_pcg(sples->Asp, sples->x, sples->b, les->rows, 250,
  158. 0.1e-9, 1);
  159. G_math_d_asum_norm(sples->x, &val, sples->rows);
  160. if ((val - (double)sples->rows) > EPSILON_ITER)
  161. {
  162. G_warning("Error in G_math_solver_sparse_pcg abs %2.20f != %i", val,
  163. sples->rows);
  164. sum++;
  165. }
  166. G_math_print_les(sples);
  167. G_math_free_les(les);
  168. G_math_free_les(sples);
  169. G_message("\t * testing pcg solver with symmetric matrix and preconditioner 2\n");
  170. les = create_normal_symmetric_les(TEST_NUM_ROWS);
  171. sples = create_sparse_symmetric_les(TEST_NUM_ROWS);
  172. G_math_solver_pcg(les->A, les->x, les->b, les->rows, 250, 0.1e-9, 2);
  173. G_math_d_asum_norm(les->x, &val, les->rows);
  174. if ((val - (double)les->rows) > EPSILON_ITER)
  175. {
  176. G_warning("Error in G_math_solver_pcg abs %2.20f != %i", val, les->rows);
  177. sum++;
  178. }
  179. G_math_print_les(les);
  180. G_math_solver_sparse_pcg(sples->Asp, sples->x, sples->b, les->rows, 250,
  181. 0.1e-9, 2);
  182. G_math_d_asum_norm(sples->x, &val, sples->rows);
  183. if ((val - (double)sples->rows) > EPSILON_ITER)
  184. {
  185. G_warning("Error in G_math_solver_sparse_pcg abs %2.20f != %i", val,
  186. sples->rows);
  187. sum++;
  188. }
  189. G_math_print_les(sples);
  190. G_math_free_les(les);
  191. G_math_free_les(sples);
  192. G_message("\t * testing pcg solver with symmetric matrix and preconditioner 3\n");
  193. les = create_normal_symmetric_les(TEST_NUM_ROWS);
  194. sples = create_sparse_symmetric_les(TEST_NUM_ROWS);
  195. G_math_solver_pcg(les->A, les->x, les->b, les->rows, 250, 0.1e-9, 3);
  196. G_math_d_asum_norm(les->x, &val, les->rows);
  197. if ((val - (double)les->rows) > EPSILON_ITER)
  198. {
  199. G_warning("Error in G_math_solver_pcg abs %2.20f != %i", val, les->rows);
  200. sum++;
  201. }
  202. G_math_print_les(les);
  203. G_math_solver_sparse_pcg(sples->Asp, sples->x, sples->b, les->rows, 250,
  204. 0.1e-9, 3);
  205. G_math_d_asum_norm(sples->x, &val, sples->rows);
  206. if ((val - (double)sples->rows) > EPSILON_ITER)
  207. {
  208. G_warning("Error in G_math_solver_sparse_pcg abs %2.20f != %i", val,
  209. sples->rows);
  210. sum++;
  211. }
  212. G_math_print_les(sples);
  213. G_math_free_les(les);
  214. G_math_free_les(sples);
  215. G_message("\t * testing cg solver with symmetric matrix\n");
  216. les = create_normal_symmetric_les(TEST_NUM_ROWS);
  217. sples = create_sparse_symmetric_les(TEST_NUM_ROWS);
  218. G_math_solver_cg(les->A, les->x, les->b, les->rows, 250, 0.1e-9);
  219. G_math_d_asum_norm(les->x, &val, les->rows);
  220. if ((val - (double)les->rows) > EPSILON_ITER)
  221. {
  222. G_warning("Error in G_math_solver_cg abs %2.20f != %i", val, les->rows);
  223. sum++;
  224. }
  225. G_math_print_les(les);
  226. G_math_solver_sparse_cg(sples->Asp, sples->x, sples->b, les->rows, 250,
  227. 0.1e-9);
  228. G_math_d_asum_norm(sples->x, &val, sples->rows);
  229. if ((val - (double)sples->rows) > EPSILON_ITER)
  230. {
  231. G_warning("Error in G_math_solver_sparse_cg abs %2.20f != %i", val,
  232. sples->rows);
  233. sum++;
  234. }
  235. G_math_print_les(sples);
  236. G_math_free_les(les);
  237. G_math_free_les(sples);
  238. G_message("\t * testing cg solver with symmetric bad conditioned matrix\n");
  239. les = create_normal_symmetric_pivot_les(TEST_NUM_ROWS);
  240. G_math_solver_cg(les->A, les->x, les->b, les->rows, 250, 0.1e-9);
  241. G_math_d_asum_norm(les->x, &val, les->rows);
  242. if ((val - (double)les->rows) > EPSILON_ITER)
  243. {
  244. G_warning("Error in G_math_solver_cg abs %2.20f != %i", val, les->rows);
  245. sum++;
  246. }
  247. G_math_print_les(les);
  248. G_math_free_les(les);
  249. G_message("\t * testing bicgstab solver with symmetric matrix\n");
  250. les = create_normal_symmetric_les(TEST_NUM_ROWS);
  251. sples = create_sparse_symmetric_les(TEST_NUM_ROWS);
  252. G_math_solver_bicgstab(les->A, les->x, les->b, les->rows, 250, 0.1e-9);
  253. G_math_d_asum_norm(les->x, &val, les->rows);
  254. if ((val - (double)les->rows) > EPSILON_ITER)
  255. {
  256. G_warning("Error in G_math_solver_bicgstab abs %2.20f != %i", val,
  257. les->rows);
  258. sum++;
  259. }
  260. G_math_print_les(les);
  261. G_math_solver_sparse_bicgstab(sples->Asp, sples->x, sples->b, les->rows,
  262. 250, 0.1e-9);
  263. G_math_d_asum_norm(sples->x, &val, sples->rows);
  264. if ((val - (double)sples->rows) > EPSILON_ITER)
  265. {
  266. G_warning("Error in G_math_solver_sparse_bicgstab abs %2.20f != %i",
  267. val, sples->rows);
  268. sum++;
  269. }
  270. G_math_print_les(sples);
  271. G_math_free_les(les);
  272. G_math_free_les(sples);
  273. G_message("\t * testing bicgstab solver with unsymmetric matrix\n");
  274. les = create_normal_unsymmetric_les(TEST_NUM_ROWS);
  275. sples = create_sparse_unsymmetric_les(TEST_NUM_ROWS);
  276. G_math_solver_bicgstab(les->A, les->x, les->b, les->rows, 250, 0.1e-9);
  277. G_math_d_asum_norm(les->x, &val, les->rows);
  278. if ((val - (double)les->rows) > EPSILON_ITER)
  279. {
  280. G_warning("Error in G_math_solver_bicgstab abs %2.20f != %i", val,
  281. les->rows);
  282. sum++;
  283. }
  284. G_math_print_les(les);
  285. G_math_solver_sparse_bicgstab(sples->Asp, sples->x, sples->b, les->rows,
  286. 250, 0.1e-9);
  287. G_math_d_asum_norm(sples->x, &val, sples->rows);
  288. if ((val - (double)les->rows) > EPSILON_ITER)
  289. {
  290. G_warning("Error in G_math_solver_sparse_bicgstab abs %2.20f != %i",
  291. val, sples->rows);
  292. sum++;
  293. }
  294. G_math_print_les(sples);
  295. G_math_free_les(les);
  296. G_math_free_les(sples);
  297. G_message("\t * testing gauss elimination solver with symmetric matrix\n");
  298. les = create_normal_symmetric_les(TEST_NUM_ROWS);
  299. G_math_solver_gauss(les->A, les->x, les->b, les->rows);
  300. G_math_d_asum_norm(les->x, &val, les->rows);
  301. if ((val - (double)les->rows) > EPSILON_DIRECT)
  302. {
  303. G_warning("Error in G_math_solver_gauss abs %2.20f != %i", val,
  304. les->rows);
  305. sum++;
  306. }
  307. G_math_print_les(les);
  308. G_math_free_les(les);
  309. G_message("\t * testing lu decomposition solver with symmetric matrix\n");
  310. les = create_normal_symmetric_les(TEST_NUM_ROWS);
  311. G_math_solver_lu(les->A, les->x, les->b, les->rows);
  312. G_math_d_asum_norm(les->x, &val, les->rows);
  313. if ((val - (double)les->rows) > EPSILON_DIRECT)
  314. {
  315. G_warning("Error in G_math_solver_gauss abs %2.20f != %i", val,
  316. les->rows);
  317. sum++;
  318. }
  319. G_math_print_les(les);
  320. G_math_free_les(les);
  321. G_message("\t * testing gauss elimination solver with unsymmetric matrix\n");
  322. les = create_normal_unsymmetric_les(TEST_NUM_ROWS);
  323. G_math_solver_gauss(les->A, les->x, les->b, les->rows);
  324. G_math_d_asum_norm(les->x, &val, les->rows);
  325. if ((val - (double)les->rows) > EPSILON_DIRECT)
  326. {
  327. G_warning("Error in G_math_solver_gauss abs %2.20f != %i", val,
  328. les->rows);
  329. sum++;
  330. }
  331. G_math_print_les(les);
  332. G_math_free_les(les);
  333. G_message("\t * testing lu decomposition solver with unsymmetric matrix\n");
  334. les = create_normal_unsymmetric_les(TEST_NUM_ROWS);
  335. G_math_solver_lu(les->A, les->x, les->b, les->rows);
  336. G_math_d_asum_norm(les->x, &val, les->rows);
  337. if ((val - (double)les->rows) > EPSILON_DIRECT)
  338. {
  339. G_warning("Error in G_math_solver_gauss abs %2.20f != %i", val,
  340. les->rows);
  341. sum++;
  342. }
  343. G_math_print_les(les);
  344. G_math_free_les(les);
  345. G_message("\t * testing gauss elimination solver with symmetric bad conditioned matrix\n");
  346. les = create_normal_symmetric_pivot_les(TEST_NUM_ROWS);
  347. G_math_solver_gauss(les->A, les->x, les->b, les->rows);
  348. G_math_d_asum_norm(les->x, &val, les->rows);
  349. if ((val - (double)les->rows) > EPSILON_DIRECT)
  350. {
  351. G_warning("Error in G_math_solver_gauss abs %2.20f != %i", val,
  352. les->rows);
  353. sum++;
  354. }
  355. G_math_print_les(les);
  356. G_math_free_les(les);
  357. G_message("\t * testing lu decomposition solver with symmetric bad conditioned matrix\n");
  358. les = create_normal_symmetric_pivot_les(TEST_NUM_ROWS);
  359. G_math_solver_lu(les->A, les->x, les->b, les->rows);
  360. G_math_d_asum_norm(les->x, &val, les->rows);
  361. if ((val - (double)les->rows) > EPSILON_DIRECT)
  362. {
  363. G_warning("Error in G_math_solver_lu abs %2.20f != %i", val,
  364. les->rows);
  365. sum++;
  366. }
  367. G_math_print_les(les);
  368. G_math_free_les(les);
  369. G_message("\t * testing cholesky decomposition solver with symmetric matrix\n");
  370. les = create_normal_symmetric_les(TEST_NUM_ROWS);
  371. /*cholesky*/G_math_solver_cholesky(les->A, les->x, les->b, les->rows,
  372. les->rows);
  373. G_math_d_asum_norm(les->x, &val, les->rows);
  374. if ((val - (double)les->rows) > EPSILON_DIRECT)
  375. {
  376. G_warning("Error in G_math_solver_solver_cholesky abs %2.20f != %i", val,
  377. les->rows);
  378. sum++;
  379. }
  380. G_math_print_les(les);
  381. G_math_free_les(les);
  382. G_message("\t * testing cholesky band decomposition solver with symmetric band matrix 1\n");
  383. les = create_normal_symmetric_les(TEST_NUM_ROWS);
  384. G_math_print_les(les);
  385. /* Create a band matrix*/
  386. G_message("\t * Creating symmetric band matrix\n");
  387. les->A = G_math_matrix_to_sband_matrix(les->A, les->rows, les->rows);
  388. G_math_print_les(les);
  389. /*cholesky*/G_math_solver_cholesky_sband(les->A, les->x, les->b, les->rows,les->rows);
  390. G_math_d_asum_norm(les->x, &val, les->rows);
  391. if ((val - (double)les->rows) > EPSILON_DIRECT)
  392. {
  393. G_warning("Error in G_math_solver_solver_cholesky_band abs %2.20f != %i", val,
  394. les->rows);
  395. sum++;
  396. }
  397. G_math_print_les(les);
  398. G_math_free_les(les);
  399. G_message("\t * testing cholesky band decomposition solver with symmetric band matrix 2\n");
  400. les = create_symmetric_band_les(TEST_NUM_ROWS);
  401. G_math_print_les(les);
  402. /*cholesky*/G_math_solver_cholesky_sband(les->A, les->x, les->b, les->rows,les->rows);
  403. G_math_d_asum_norm(les->x, &val, les->rows);
  404. if ((val - (double)les->rows) > EPSILON_DIRECT)
  405. {
  406. G_warning("Error in G_math_solver_solver_cholesky_band abs %2.20f != %i", val,
  407. les->rows);
  408. sum++;
  409. }
  410. G_math_print_les(les);
  411. G_math_free_les(les);
  412. G_message("\t * testing cg solver with symmetric band matrix\n");
  413. les = create_symmetric_band_les(TEST_NUM_ROWS);
  414. G_math_solver_cg_sband(les->A, les->x, les->b, les->rows, les->rows, 250, 0.1e-9);
  415. G_math_d_asum_norm(les->x, &val, les->rows);
  416. if ((val - (double)les->rows) > EPSILON_ITER)
  417. {
  418. G_warning("Error in G_math_solver_cg_sband abs %2.20f != %i", val, les->rows);
  419. sum++;
  420. }
  421. G_math_print_les(les);
  422. G_math_free_les(les);
  423. return sum;
  424. }