test_ccmath_wrapper.c 6.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231
  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_ccmath_wrapper(void);
  26. /* ************************************************************************* */
  27. /* Performe the solver unit tests ****************************************** */
  28. /* ************************************************************************* */
  29. int unit_test_ccmath_wrapper(void)
  30. {
  31. int sum = 0;
  32. G_message(_("\n++ Running ccmath wrapper unit tests ++"));
  33. sum += test_ccmath_wrapper();
  34. if (sum > 0)
  35. G_warning(_("\n-- ccmath wrapper unit tests failure --"));
  36. else
  37. G_message(_("\n-- ccmath wrapper unit tests finished successfully --"));
  38. return sum;
  39. }
  40. /* *************************************************************** */
  41. /* Test all implemented ccmath wrapper *** */
  42. /* *************************************************************** */
  43. int test_ccmath_wrapper(void)
  44. {
  45. G_math_les *les;
  46. int sum = 0;
  47. double val = 0.0, val2 = 0.0;
  48. G_message("\t * testing ccmath lu solver with symmetric matrix\n");
  49. les = create_normal_symmetric_les(TEST_NUM_ROWS);
  50. G_math_d_copy(les->b, les->x, les->rows);
  51. G_math_solv(les->A, les->x, les->rows);
  52. G_math_print_les(les);
  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_solv abs %2.20f != %i", val,
  57. les->rows);
  58. sum++;
  59. }
  60. G_math_free_les(les);
  61. G_message("\t * testing ccmath lu solver with unsymmetric matrix\n");
  62. les = create_normal_unsymmetric_les(TEST_NUM_ROWS);
  63. G_math_d_copy(les->b, les->x, les->rows);
  64. G_math_solvps(les->A, les->x, les->rows);
  65. G_math_print_les(les);
  66. G_math_d_asum_norm(les->x, &val, les->rows);
  67. if ((val - (double)les->rows) > EPSILON_ITER)
  68. {
  69. G_warning("Error in G_math_solv abs %2.20f != %i", val,
  70. les->rows);
  71. sum++;
  72. }
  73. G_math_free_les(les);
  74. G_message("\t * testing ccmath positive definite solver with symmetric matrix\n");
  75. les = create_normal_symmetric_les(TEST_NUM_ROWS);
  76. G_math_d_copy(les->b, les->x, les->rows);
  77. G_math_solvps(les->A, les->x, les->rows);
  78. G_math_print_les(les);
  79. G_math_d_asum_norm(les->x, &val, les->rows);
  80. if ((val - (double)les->rows) > EPSILON_ITER)
  81. {
  82. G_warning("Error in G_math_solvps abs %2.20f != %i", val,
  83. les->rows);
  84. sum++;
  85. }
  86. G_math_free_les(les);
  87. G_message("\t * testing ccmath matrix inversion with symmetric matrix\n");
  88. les = create_normal_symmetric_les(TEST_NUM_ROWS);
  89. G_math_minv(les->A, les->rows);
  90. G_math_d_Ax(les->A, les->b, les->x, les->rows, les->rows);
  91. G_math_print_les(les);
  92. G_math_d_asum_norm(les->x, &val, les->rows);
  93. if ((val - (double)les->rows) > EPSILON_ITER)
  94. {
  95. G_warning("Error in G_math_minv abs %2.20f != %i", val,
  96. les->rows);
  97. sum++;
  98. }
  99. G_math_free_les(les);
  100. G_message("\t * testing ccmath matrix inversion with unsymmetric matrix\n");
  101. les = create_normal_unsymmetric_les(TEST_NUM_ROWS);
  102. G_math_minv(les->A, les->rows);
  103. G_math_d_Ax(les->A, les->b, les->x, les->rows, les->rows);
  104. G_math_print_les(les);
  105. G_math_d_asum_norm(les->x, &val, les->rows);
  106. if ((val - (double)les->rows) > EPSILON_ITER)
  107. {
  108. G_warning("Error in G_math_minv abs %2.20f != %i", val,
  109. les->rows);
  110. sum++;
  111. }
  112. G_math_free_les(les);
  113. G_message("\t * testing ccmath positive definite matrix inversion with symmetric matrix\n");
  114. les = create_normal_symmetric_les(TEST_NUM_ROWS);
  115. G_math_psinv(les->A, les->rows);
  116. G_math_d_Ax(les->A, les->b, les->x, les->rows, les->rows);
  117. G_math_print_les(les);
  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_psinv abs %2.20f != %i", val,
  122. les->rows);
  123. sum++;
  124. }
  125. G_math_free_les(les);
  126. G_message("\t * testing ccmath eigenvalue solver with symmetric matrix\n");
  127. les = create_normal_symmetric_les(TEST_NUM_ROWS);
  128. // Results of the eigenvalue computation with ocatve
  129. les->b[9] = 0.043264;
  130. les->b[8] = 0.049529;
  131. les->b[7] = 0.057406;
  132. les->b[6] = 0.067696;
  133. les->b[5] = 0.081639;
  134. les->b[4] = 0.101357;
  135. les->b[3] = 0.130298;
  136. les->b[2] = 0.174596;
  137. les->b[1] = 0.256157;
  138. les->b[0] = 0.502549;
  139. G_math_eigval(les->A, les->x, les->rows);
  140. G_math_print_les(les);
  141. G_math_d_asum_norm(les->b, &val, les->rows);
  142. G_math_d_asum_norm(les->x, &val2, les->rows);
  143. if ((val - val2) > EPSILON_ITER)
  144. {
  145. G_warning("Error in G_math_eigv abs %2.20f != %f", val,
  146. val2);
  147. sum++;
  148. }
  149. G_math_free_les(les);
  150. G_message("\t * testing ccmath eigenvector computation with symmetric matrix\n");
  151. les = create_normal_symmetric_les(TEST_NUM_ROWS);
  152. // Results of the eigenvalue computation with ocatve
  153. les->b[9] = 0.043264;
  154. les->b[8] = 0.049529;
  155. les->b[7] = 0.057406;
  156. les->b[6] = 0.067696;
  157. les->b[5] = 0.081639;
  158. les->b[4] = 0.101357;
  159. les->b[3] = 0.130298;
  160. les->b[2] = 0.174596;
  161. les->b[1] = 0.256157;
  162. les->b[0] = 0.502549;
  163. G_math_eigen(les->A, les->x, les->rows);
  164. G_math_print_les(les);
  165. G_math_d_asum_norm(les->b, &val, les->rows);
  166. G_math_d_asum_norm(les->x, &val2, les->rows);
  167. if ((val - val2) > EPSILON_ITER)
  168. {
  169. G_warning("Error in G_math_eigen abs %2.20f != %f", val,
  170. val2);
  171. sum++;
  172. }
  173. G_math_free_les(les);
  174. G_message("\t * testing ccmath singulare value decomposition with symmetric matrix\n");
  175. les = create_normal_symmetric_les(TEST_NUM_ROWS);
  176. G_math_svdval(les->x, les->A, les->rows, les->rows);
  177. G_math_print_les(les);
  178. G_math_free_les(les);
  179. return sum;
  180. }