bench_solver_krylov.c 4.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123
  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: benchmarking the krylov subspace solvers
  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. #include <sys/time.h>
  23. /* prototypes */
  24. static int bench_solvers(int rows);
  25. /* ************************************************************************* */
  26. /* Performe the solver unit tests ****************************************** */
  27. /* ************************************************************************* */
  28. int bench_solvers_krylov(int rows) {
  29. G_message(_("\n++ Running krylov solver benchmark ++"));
  30. bench_solvers(rows);
  31. return 1;
  32. }
  33. /* *************************************************************** */
  34. /* Test all implemented solvers for sparse and normal matrix *** */
  35. /* *************************************************************** */
  36. int bench_solvers(int rows) {
  37. G_math_les *les;
  38. G_math_les *sples;
  39. struct timeval tstart;
  40. struct timeval tend;
  41. G_message("\t * benchmarking pcg solver with symmetric matrix and preconditioner 1\n");
  42. les = create_normal_symmetric_les(rows);
  43. sples = create_sparse_symmetric_les(rows);
  44. gettimeofday(&tstart, NULL);
  45. G_math_solver_pcg(les->A, les->x, les->b, les->rows, 250, 0.1e-9, 1);
  46. gettimeofday(&tend, NULL);
  47. printf("Computation time pcg normal matrix: %g\n", compute_time_difference(tstart, tend));
  48. gettimeofday(&tstart, NULL);
  49. G_math_solver_sparse_pcg(sples->Asp, sples->x, sples->b, les->rows, 250,
  50. 0.1e-9, 1);
  51. gettimeofday(&tend, NULL);
  52. printf("Computation time pcg sparse matrix: %g\n", compute_time_difference(tstart, tend));
  53. G_math_free_les(les);
  54. G_math_free_les(sples);
  55. G_message("\t * benchmark cg solver with symmetric matrix\n");
  56. les = create_normal_symmetric_les(rows);
  57. sples = create_sparse_symmetric_les(rows);
  58. gettimeofday(&tstart, NULL);
  59. G_math_solver_cg(les->A, les->x, les->b, les->rows, 250, 0.1e-9);
  60. gettimeofday(&tend, NULL);
  61. printf("Computation time cg normal matrix: %g\n", compute_time_difference(tstart, tend));
  62. gettimeofday(&tstart, NULL);
  63. G_math_solver_sparse_cg(sples->Asp, sples->x, sples->b, les->rows, 250,
  64. 0.1e-9);
  65. gettimeofday(&tend, NULL);
  66. printf("Computation time cg sparse matrix: %g\n", compute_time_difference(tstart, tend));
  67. G_math_free_les(les);
  68. G_math_free_les(sples);
  69. G_message("\t * benchmark cg solver with symmetric band matrix\n");
  70. les = create_symmetric_band_les(rows);
  71. gettimeofday(&tstart, NULL);
  72. G_math_solver_cg_sband(les->A, les->x, les->b, les->rows, les->rows, 250, 0.1e-9);
  73. gettimeofday(&tend, NULL);
  74. printf("Computation time cg symmetric band matrix: %g\n", compute_time_difference(tstart, tend));
  75. G_math_free_les(les);
  76. G_message("\t * benchmark bicgstab solver with unsymmetric matrix\n");
  77. les = create_normal_unsymmetric_les(rows);
  78. sples = create_sparse_unsymmetric_les(rows);
  79. gettimeofday(&tstart, NULL);
  80. G_math_solver_bicgstab(les->A, les->x, les->b, les->rows, 250, 0.1e-9);
  81. gettimeofday(&tend, NULL);
  82. printf("Computation time bicgstab normal matrix: %g\n", compute_time_difference(tstart, tend));
  83. gettimeofday(&tstart, NULL);
  84. G_math_solver_sparse_bicgstab(sples->Asp, sples->x, sples->b, les->rows,
  85. 250, 0.1e-9);
  86. gettimeofday(&tend, NULL);
  87. printf("Computation time bicgstab sparse matrix: %g\n", compute_time_difference(tstart, tend));
  88. G_math_free_les(les);
  89. G_math_free_les(sples);
  90. return 1;
  91. }