test_gwflow.c 17 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534
  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: gwflow integration tests
  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 <grass/gmath.h>
  17. #include <grass/gis.h>
  18. #include <grass/N_pde.h>
  19. #include <grass/gmath.h>
  20. #include <grass/N_gwflow.h>
  21. #include "test_gpde_lib.h"
  22. /*redefine */
  23. #define TEST_N_NUM_DEPTHS_LOCAL 2
  24. #define TEST_N_NUM_ROWS_LOCAL 3
  25. #define TEST_N_NUM_COLS_LOCAL 3
  26. /* prototypes */
  27. static N_gwflow_data2d *create_gwflow_data_2d(void);
  28. static N_gwflow_data3d *create_gwflow_data_3d(void);
  29. static int test_gwflow_2d(void);
  30. static int test_gwflow_3d(void);
  31. /* *************************************************************** */
  32. /* Performe the gwflow integration tests ************************* */
  33. /* *************************************************************** */
  34. int integration_test_gwflow(void)
  35. {
  36. int sum = 0;
  37. G_message("\n++ Running gwflow integration tests ++");
  38. G_message("\t 1. testing 2d gwflow");
  39. sum += test_gwflow_2d();
  40. G_message("\t 2. testing 3d gwflow");
  41. sum += test_gwflow_3d();
  42. if (sum > 0)
  43. G_warning("\n-- gwflow integration tests failure --");
  44. else
  45. G_message("\n-- gwflow integration tests finished successfully --");
  46. return sum;
  47. }
  48. /* *************************************************************** */
  49. /* Create valid groundwater flow data **************************** */
  50. /* *************************************************************** */
  51. N_gwflow_data3d *create_gwflow_data_3d(void)
  52. {
  53. N_gwflow_data3d *data;
  54. int i, j, k;
  55. data =
  56. N_alloc_gwflow_data3d(TEST_N_NUM_COLS_LOCAL, TEST_N_NUM_ROWS_LOCAL,
  57. TEST_N_NUM_DEPTHS_LOCAL, 1, 1);
  58. #pragma omp parallel for private (i, j, k) shared (data)
  59. for (k = 0; k < TEST_N_NUM_DEPTHS_LOCAL; k++)
  60. for (j = 0; j < TEST_N_NUM_ROWS_LOCAL; j++) {
  61. for (i = 0; i < TEST_N_NUM_COLS_LOCAL; i++) {
  62. if (j == 0) {
  63. N_put_array_3d_d_value(data->phead, i, j, k, 50);
  64. N_put_array_3d_d_value(data->phead_start, i, j, k, 50);
  65. N_put_array_3d_d_value(data->status, i, j, k, 2);
  66. }
  67. else {
  68. N_put_array_3d_d_value(data->phead, i, j, k, 40);
  69. N_put_array_3d_d_value(data->phead_start, i, j, k, 40);
  70. N_put_array_3d_d_value(data->status, i, j, k, 1);
  71. }
  72. N_put_array_3d_d_value(data->hc_x, i, j, k, 0.0001);
  73. N_put_array_3d_d_value(data->hc_y, i, j, k, 0.0001);
  74. N_put_array_3d_d_value(data->hc_z, i, j, k, 0.0001);
  75. N_put_array_3d_d_value(data->q, i, j, k, 0.0);
  76. N_put_array_3d_d_value(data->s, i, j, k, 0.001);
  77. N_put_array_2d_d_value(data->r, i, j, 0.0);
  78. N_put_array_3d_d_value(data->nf, i, j, k, 0.1);
  79. }
  80. }
  81. return data;
  82. }
  83. /* *************************************************************** */
  84. /* Create valid groundwater flow data **************************** */
  85. /* *************************************************************** */
  86. N_gwflow_data2d *create_gwflow_data_2d(void)
  87. {
  88. int i, j;
  89. N_gwflow_data2d *data;
  90. data = N_alloc_gwflow_data2d(TEST_N_NUM_COLS_LOCAL, TEST_N_NUM_ROWS_LOCAL, 1, 1);
  91. #pragma omp parallel for private (i, j) shared (data)
  92. for (j = 0; j < TEST_N_NUM_ROWS_LOCAL; j++) {
  93. for (i = 0; i < TEST_N_NUM_COLS_LOCAL; i++) {
  94. if (j == 0) {
  95. N_put_array_2d_d_value(data->phead, i, j, 50);
  96. N_put_array_2d_d_value(data->phead_start, i, j, 50);
  97. N_put_array_2d_d_value(data->status, i, j, 2);
  98. }
  99. else {
  100. N_put_array_2d_d_value(data->phead, i, j, 40);
  101. N_put_array_2d_d_value(data->phead_start, i, j, 40);
  102. N_put_array_2d_d_value(data->status, i, j, 1);
  103. }
  104. N_put_array_2d_d_value(data->hc_x, i, j, 30.0001);
  105. N_put_array_2d_d_value(data->hc_y, i, j, 30.0001);
  106. N_put_array_2d_d_value(data->q, i, j, 0.0);
  107. N_put_array_2d_d_value(data->s, i, j, 0.001);
  108. N_put_array_2d_d_value(data->r, i, j, 0.0);
  109. N_put_array_2d_d_value(data->nf, i, j, 0.1);
  110. N_put_array_2d_d_value(data->top, i, j, 20.0);
  111. N_put_array_2d_d_value(data->bottom, i, j, 0.0);
  112. }
  113. }
  114. return data;
  115. }
  116. /* *************************************************************** */
  117. /* Test the groundwater flow in 3d with different solvers ******** */
  118. /* *************************************************************** */
  119. int test_gwflow_3d(void)
  120. {
  121. N_gwflow_data3d *data;
  122. N_geom_data *geom;
  123. N_les *les;
  124. N_les_callback_3d *call;
  125. call = N_alloc_les_callback_3d();
  126. N_set_les_callback_3d_func(call, (*N_callback_gwflow_3d)); /*gwflow 3d */
  127. data = create_gwflow_data_3d();
  128. data->dt = 86400;
  129. geom = N_alloc_geom_data();
  130. geom->dx = 10;
  131. geom->dy = 15;
  132. geom->dz = 3;
  133. geom->Az = 150;
  134. geom->depths = TEST_N_NUM_DEPTHS_LOCAL;
  135. geom->rows = TEST_N_NUM_ROWS_LOCAL;
  136. geom->cols = TEST_N_NUM_COLS_LOCAL;
  137. /*Assemble the matrix */
  138. /*
  139. */
  140. /*CG*/ les =
  141. N_assemble_les_3d(N_SPARSE_LES, geom, data->status, data->phead_start,
  142. (void *)data, call);
  143. G_math_solver_sparse_cg(les->Asp, les->x, les->b, les->rows, 100, 0.1e-8);
  144. N_print_les(les);
  145. N_free_les(les);
  146. /*PCG G_MATH_DIAGONAL_PRECONDITION*/ les =
  147. N_assemble_les_3d(N_SPARSE_LES, geom, data->status, data->phead_start,
  148. (void *)data, call);
  149. G_math_solver_sparse_pcg(les->Asp, les->x, les->b, les->rows, 100, 0.1e-8, G_MATH_DIAGONAL_PRECONDITION);
  150. N_print_les(les);
  151. N_free_les(les);
  152. /*PCG G_MATH_ROWSCALE_EUKLIDNORM_PRECONDITION*/ les =
  153. N_assemble_les_3d(N_SPARSE_LES, geom, data->status, data->phead_start,
  154. (void *)data, call);
  155. G_math_solver_sparse_pcg(les->Asp, les->x, les->b, les->rows, 100, 0.1e-8, G_MATH_ROWSCALE_EUKLIDNORM_PRECONDITION);
  156. N_print_les(les);
  157. N_free_les(les);
  158. /*PCG G_MATH_ROWSCALE_ABSSUMNORM_PRECONDITION*/ les =
  159. N_assemble_les_3d(N_SPARSE_LES, geom, data->status, data->phead_start,
  160. (void *)data, call);
  161. G_math_solver_sparse_pcg(les->Asp, les->x, les->b, les->rows, 100, 0.1e-8, G_MATH_ROWSCALE_ABSSUMNORM_PRECONDITION);
  162. N_print_les(les);
  163. N_free_les(les);
  164. /*CG*/ les =
  165. N_assemble_les_3d_dirichlet(N_SPARSE_LES, geom, data->status, data->phead_start,
  166. (void *)data, call);
  167. N_les_integrate_dirichlet_3d(les, geom, data->status, data->phead_start);
  168. G_math_solver_sparse_cg(les->Asp, les->x, les->b, les->rows, 100, 0.1e-8);
  169. N_print_les(les);
  170. N_free_les(les);
  171. /*CG*/ les =
  172. N_assemble_les_3d(N_NORMAL_LES, geom, data->status, data->phead_start,
  173. (void *)data, call);
  174. G_math_solver_cg(les->A, les->x, les->b, les->rows, 100, 0.1e-8);
  175. N_print_les(les);
  176. N_free_les(les);
  177. /*PCG G_MATH_DIAGONAL_PRECONDITION*/ les =
  178. N_assemble_les_3d(N_NORMAL_LES, geom, data->status, data->phead_start,
  179. (void *)data, call);
  180. G_math_solver_pcg(les->A, les->x, les->b, les->rows, 100, 0.1e-8, G_MATH_DIAGONAL_PRECONDITION);
  181. N_print_les(les);
  182. N_free_les(les);
  183. /*PCG G_MATH_ROWSCALE_EUKLIDNORM_PRECONDITION*/ les =
  184. N_assemble_les_3d(N_NORMAL_LES, geom, data->status, data->phead_start,
  185. (void *)data, call);
  186. G_math_solver_pcg(les->A, les->x, les->b, les->rows, 100, 0.1e-8, G_MATH_ROWSCALE_EUKLIDNORM_PRECONDITION);
  187. N_print_les(les);
  188. N_free_les(les);
  189. /*PCG G_MATH_ROWSCALE_ABSSUMNORM_PRECONDITION*/ les =
  190. N_assemble_les_3d(N_NORMAL_LES, geom, data->status, data->phead_start,
  191. (void *)data, call);
  192. G_math_solver_pcg(les->A, les->x, les->b, les->rows, 100, 0.1e-8, G_MATH_ROWSCALE_ABSSUMNORM_PRECONDITION);
  193. N_print_les(les);
  194. N_free_les(les);
  195. /*CG*/ les =
  196. N_assemble_les_3d_dirichlet(N_NORMAL_LES, geom, data->status, data->phead_start,
  197. (void *)data, call);
  198. N_les_integrate_dirichlet_3d(les, geom, data->status, data->phead_start);
  199. G_math_solver_cg(les->A, les->x, les->b, les->rows, 100, 0.1e-8);
  200. N_print_les(les);
  201. N_free_les(les);
  202. /*BICG*/ les =
  203. N_assemble_les_3d(N_SPARSE_LES, geom, data->status, data->phead_start,
  204. (void *)data, call);
  205. G_math_solver_sparse_bicgstab(les->Asp, les->x, les->b, les->rows, 100, 0.1e-8);
  206. N_print_les(les);
  207. N_free_les(les);
  208. /*BICG*/ les =
  209. N_assemble_les_3d(N_NORMAL_LES, geom, data->status, data->phead_start,
  210. (void *)data, call);
  211. G_math_solver_bicgstab(les->A, les->x, les->b, les->rows, 100, 0.1e-8);
  212. N_print_les(les);
  213. N_free_les(les);
  214. /*BICG*/ les =
  215. N_assemble_les_3d_dirichlet(N_SPARSE_LES, geom, data->status, data->phead_start,
  216. (void *)data, call);
  217. N_les_integrate_dirichlet_3d(les, geom, data->status, data->phead_start);
  218. G_math_solver_sparse_bicgstab(les->Asp, les->x, les->b, les->rows, 100, 0.1e-8);
  219. N_print_les(les);
  220. N_free_les(les);
  221. /*BICG*/ les =
  222. N_assemble_les_3d_dirichlet(N_NORMAL_LES, geom, data->status, data->phead_start,
  223. (void *)data, call);
  224. N_les_integrate_dirichlet_3d(les, geom, data->status, data->phead_start);
  225. G_math_solver_bicgstab(les->A, les->x, les->b, les->rows, 100, 0.1e-8);
  226. N_print_les(les);
  227. N_free_les(les);
  228. /*GUASS*/ les =
  229. N_assemble_les_3d(N_NORMAL_LES, geom, data->status, data->phead_start,
  230. (void *)data, call);
  231. G_math_solver_gauss(les->A, les->x, les->b, les->rows);
  232. N_print_les(les);
  233. N_free_les(les);
  234. /*LU*/ les =
  235. N_assemble_les_3d(N_NORMAL_LES, geom, data->status, data->phead_start,
  236. (void *)data, call);
  237. G_math_solver_lu(les->A, les->x, les->b, les->rows);
  238. N_print_les(les);
  239. N_free_les(les);
  240. /*GUASS*/ les =
  241. N_assemble_les_3d_dirichlet(N_NORMAL_LES, geom, data->status, data->phead_start,
  242. (void *)data, call);
  243. N_les_integrate_dirichlet_3d(les, geom, data->status, data->phead_start);
  244. G_math_solver_gauss(les->A, les->x, les->b, les->rows);
  245. N_print_les(les);
  246. N_free_les(les);
  247. /*LU*/ les =
  248. N_assemble_les_3d_dirichlet(N_NORMAL_LES, geom, data->status, data->phead_start,
  249. (void *)data, call);
  250. N_les_integrate_dirichlet_3d(les, geom, data->status, data->phead_start);
  251. G_math_solver_lu(les->A, les->x, les->b, les->rows);
  252. N_print_les(les);
  253. N_free_les(les);
  254. /*Cholesky*/ les =
  255. N_assemble_les_3d(N_NORMAL_LES, geom, data->status, data->phead_start,
  256. (void *)data, call);
  257. G_math_solver_cholesky(les->A, les->x, les->b, les->rows, les->rows);
  258. N_print_les(les);
  259. N_free_les(les);
  260. /*Cholesky*/ les =
  261. N_assemble_les_3d_dirichlet(N_NORMAL_LES, geom, data->status, data->phead_start,
  262. (void *)data, call);
  263. N_les_integrate_dirichlet_3d(les, geom, data->status, data->phead_start);
  264. G_math_solver_cholesky(les->A, les->x, les->b, les->rows, les->rows);
  265. N_print_les(les);
  266. N_free_les(les);
  267. N_free_gwflow_data3d(data);
  268. G_free(geom);
  269. G_free(call);
  270. return 0;
  271. }
  272. /* *************************************************************** */
  273. int test_gwflow_2d(void)
  274. {
  275. N_gwflow_data2d *data;
  276. N_geom_data *geom;
  277. N_les *les;
  278. N_les_callback_2d *call;
  279. /*set the callback */
  280. call = N_alloc_les_callback_2d();
  281. N_set_les_callback_2d_func(call, (*N_callback_gwflow_2d));
  282. data = create_gwflow_data_2d();
  283. data->dt = 600;
  284. geom = N_alloc_geom_data();
  285. geom->dx = 10;
  286. geom->dy = 15;
  287. geom->Az = 150;
  288. geom->rows = TEST_N_NUM_ROWS_LOCAL;
  289. geom->cols = TEST_N_NUM_COLS_LOCAL;
  290. /*Assemble the matrix */
  291. /*
  292. */
  293. /*CG*/ les =
  294. N_assemble_les_2d(N_SPARSE_LES, geom, data->status, data->phead_start,
  295. (void *)data, call);
  296. G_math_solver_sparse_cg(les->Asp, les->x, les->b, les->rows, 100, 0.1e-8);
  297. N_print_les(les);
  298. N_free_les(les);
  299. /*PCG G_MATH_DIAGONAL_PRECONDITION*/ les =
  300. N_assemble_les_2d(N_SPARSE_LES, geom, data->status, data->phead_start,
  301. (void *)data, call);
  302. G_math_solver_sparse_pcg(les->Asp, les->x, les->b, les->rows, 100, 0.1e-8, G_MATH_DIAGONAL_PRECONDITION);
  303. N_print_les(les);
  304. N_free_les(les);
  305. /*PCG G_MATH_ROWSCALE_EUKLIDNORM_PRECONDITION*/ les =
  306. N_assemble_les_2d(N_SPARSE_LES, geom, data->status, data->phead_start,
  307. (void *)data, call);
  308. G_math_solver_sparse_pcg(les->Asp, les->x, les->b, les->rows, 100, 0.1e-8, G_MATH_ROWSCALE_EUKLIDNORM_PRECONDITION);
  309. N_print_les(les);
  310. N_free_les(les);
  311. /*PCG G_MATH_ROWSCALE_ABSSUMNORM_PRECONDITION*/ les =
  312. N_assemble_les_2d(N_SPARSE_LES, geom, data->status, data->phead_start,
  313. (void *)data, call);
  314. G_math_solver_sparse_pcg(les->Asp, les->x, les->b, les->rows, 100, 0.1e-8, G_MATH_ROWSCALE_ABSSUMNORM_PRECONDITION);
  315. N_print_les(les);
  316. N_free_les(les);
  317. /*CG*/ les =
  318. N_assemble_les_2d_dirichlet(N_SPARSE_LES, geom, data->status, data->phead_start,
  319. (void *)data, call);
  320. N_les_integrate_dirichlet_2d(les, geom, data->status, data->phead_start);
  321. G_math_solver_sparse_cg(les->Asp, les->x, les->b, les->rows, 100, 0.1e-8);
  322. N_print_les(les);
  323. N_free_les(les);
  324. /*CG*/ les =
  325. N_assemble_les_2d(N_NORMAL_LES, geom, data->status, data->phead_start,
  326. (void *)data, call);
  327. G_math_solver_cg(les->A, les->x, les->b, les->rows, 100, 0.1e-8);
  328. N_print_les(les);
  329. N_free_les(les);
  330. /*PCG G_MATH_DIAGONAL_PRECONDITION*/ les =
  331. N_assemble_les_2d(N_NORMAL_LES, geom, data->status, data->phead_start,
  332. (void *)data, call);
  333. G_math_solver_pcg(les->A, les->x, les->b, les->rows, 100, 0.1e-8, G_MATH_DIAGONAL_PRECONDITION);
  334. N_print_les(les);
  335. N_free_les(les);
  336. /*PCG G_MATH_ROWSCALE_EUKLIDNORM_PRECONDITION*/ les =
  337. N_assemble_les_2d(N_NORMAL_LES, geom, data->status, data->phead_start,
  338. (void *)data, call);
  339. G_math_solver_pcg(les->A, les->x, les->b, les->rows, 100, 0.1e-8, G_MATH_ROWSCALE_EUKLIDNORM_PRECONDITION);
  340. N_print_les(les);
  341. N_free_les(les);
  342. /*PCG G_MATH_ROWSCALE_ABSSUMNORM_PRECONDITION*/ les =
  343. N_assemble_les_2d(N_NORMAL_LES, geom, data->status, data->phead_start,
  344. (void *)data, call);
  345. G_math_solver_pcg(les->A, les->x, les->b, les->rows, 100, 0.1e-8, G_MATH_ROWSCALE_ABSSUMNORM_PRECONDITION);
  346. N_print_les(les);
  347. N_free_les(les);
  348. /*CG*/ les =
  349. N_assemble_les_2d_dirichlet(N_NORMAL_LES, geom, data->status, data->phead_start,
  350. (void *)data, call);
  351. N_les_integrate_dirichlet_2d(les, geom, data->status, data->phead_start);
  352. G_math_solver_cg(les->A, les->x, les->b, les->rows, 100, 0.1e-8);
  353. N_print_les(les);
  354. N_free_les(les);
  355. /*BICG*/ les =
  356. N_assemble_les_2d(N_SPARSE_LES, geom, data->status, data->phead_start,
  357. (void *)data, call);
  358. G_math_solver_sparse_bicgstab(les->Asp, les->x, les->b, les->rows, 100, 0.1e-8);
  359. N_print_les(les);
  360. N_free_les(les);
  361. /*BICG*/ les =
  362. N_assemble_les_2d(N_NORMAL_LES, geom, data->status, data->phead_start,
  363. (void *)data, call);
  364. G_math_solver_bicgstab(les->A, les->x, les->b, les->rows, 100, 0.1e-8);
  365. N_print_les(les);
  366. N_free_les(les);
  367. /*BICG*/ les =
  368. N_assemble_les_2d_dirichlet(N_SPARSE_LES, geom, data->status, data->phead_start,
  369. (void *)data, call);
  370. N_les_integrate_dirichlet_2d(les, geom, data->status, data->phead_start);
  371. G_math_solver_sparse_bicgstab(les->Asp, les->x, les->b, les->rows, 100, 0.1e-8);
  372. N_print_les(les);
  373. N_free_les(les);
  374. /*BICG*/ les =
  375. N_assemble_les_2d_dirichlet(N_NORMAL_LES, geom, data->status, data->phead_start,
  376. (void *)data, call);
  377. N_les_integrate_dirichlet_2d(les, geom, data->status, data->phead_start);
  378. G_math_solver_bicgstab(les->A, les->x, les->b, les->rows, 100, 0.1e-8);
  379. N_print_les(les);
  380. N_free_les(les);
  381. /*GUASS*/ les =
  382. N_assemble_les_2d(N_NORMAL_LES, geom, data->status, data->phead_start,
  383. (void *)data, call);
  384. G_math_solver_gauss(les->A, les->x, les->b, les->rows);
  385. N_print_les(les);
  386. N_free_les(les);
  387. /*LU*/ les =
  388. N_assemble_les_2d(N_NORMAL_LES, geom, data->status, data->phead_start,
  389. (void *)data, call);
  390. G_math_solver_lu(les->A, les->x, les->b, les->rows);
  391. N_print_les(les);
  392. N_free_les(les);
  393. /*GUASS*/ les =
  394. N_assemble_les_2d_dirichlet(N_NORMAL_LES, geom, data->status, data->phead_start,
  395. (void *)data, call);
  396. N_les_integrate_dirichlet_2d(les, geom, data->status, data->phead_start);
  397. G_math_solver_gauss(les->A, les->x, les->b, les->rows);
  398. N_print_les(les);
  399. N_free_les(les);
  400. /*LU*/ les =
  401. N_assemble_les_2d_dirichlet(N_NORMAL_LES, geom, data->status, data->phead_start,
  402. (void *)data, call);
  403. N_les_integrate_dirichlet_2d(les, geom, data->status, data->phead_start);
  404. G_math_solver_lu(les->A, les->x, les->b, les->rows);
  405. N_print_les(les);
  406. N_free_les(les);
  407. /*Cholesky*/ les =
  408. N_assemble_les_2d(N_NORMAL_LES, geom, data->status, data->phead_start,
  409. (void *)data, call);
  410. G_math_solver_cholesky(les->A, les->x, les->b, les->rows, les->rows);
  411. N_print_les(les);
  412. N_free_les(les);
  413. /*Cholesky*/ les =
  414. N_assemble_les_2d_dirichlet(N_NORMAL_LES, geom, data->status, data->phead_start,
  415. (void *)data, call);
  416. N_les_integrate_dirichlet_2d(les, geom, data->status, data->phead_start);
  417. G_math_solver_cholesky(les->A, les->x, les->b, les->rows, les->rows);
  418. N_print_les(les);
  419. N_free_les(les);
  420. N_free_gwflow_data2d(data);
  421. G_free(geom);
  422. G_free(call);
  423. return 0;
  424. }