test_gwflow.c 17 KB

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