diff options
Diffstat (limited to 'lib/Numeric/GSL/gsl-aux.c')
-rw-r--r-- | lib/Numeric/GSL/gsl-aux.c | 62 |
1 files changed, 32 insertions, 30 deletions
diff --git a/lib/Numeric/GSL/gsl-aux.c b/lib/Numeric/GSL/gsl-aux.c index 5c17836..7f1cf68 100644 --- a/lib/Numeric/GSL/gsl-aux.c +++ b/lib/Numeric/GSL/gsl-aux.c | |||
@@ -32,7 +32,7 @@ | |||
32 | #include <gsl/gsl_complex_math.h> | 32 | #include <gsl/gsl_complex_math.h> |
33 | #include <gsl/gsl_rng.h> | 33 | #include <gsl/gsl_rng.h> |
34 | #include <gsl/gsl_randist.h> | 34 | #include <gsl/gsl_randist.h> |
35 | #include <gsl/gsl_odeiv.h> | 35 | #include <gsl/gsl_odeiv2.h> |
36 | #include <gsl/gsl_multifit_nlin.h> | 36 | #include <gsl/gsl_multifit_nlin.h> |
37 | #include <string.h> | 37 | #include <string.h> |
38 | #include <stdio.h> | 38 | #include <stdio.h> |
@@ -1356,38 +1356,38 @@ int ode(int method, double h, double eps_abs, double eps_rel, | |||
1356 | int jac(double, int, const double*, int, int, double*), | 1356 | int jac(double, int, const double*, int, int, double*), |
1357 | KRVEC(xi), KRVEC(ts), RMAT(sol)) { | 1357 | KRVEC(xi), KRVEC(ts), RMAT(sol)) { |
1358 | 1358 | ||
1359 | const gsl_odeiv_step_type * T; | 1359 | const gsl_odeiv2_step_type * T; |
1360 | 1360 | ||
1361 | switch(method) { | 1361 | switch(method) { |
1362 | case 0 : {T = gsl_odeiv_step_rk2; break; } | 1362 | case 0 : {T = gsl_odeiv2_step_rk2; break; } |
1363 | case 1 : {T = gsl_odeiv_step_rk4; break; } | 1363 | case 1 : {T = gsl_odeiv2_step_rk4; break; } |
1364 | case 2 : {T = gsl_odeiv_step_rkf45; break; } | 1364 | case 2 : {T = gsl_odeiv2_step_rkf45; break; } |
1365 | case 3 : {T = gsl_odeiv_step_rkck; break; } | 1365 | case 3 : {T = gsl_odeiv2_step_rkck; break; } |
1366 | case 4 : {T = gsl_odeiv_step_rk8pd; break; } | 1366 | case 4 : {T = gsl_odeiv2_step_rk8pd; break; } |
1367 | case 5 : {T = gsl_odeiv_step_rk2imp; break; } | 1367 | case 5 : {T = gsl_odeiv2_step_rk2imp; break; } |
1368 | case 6 : {T = gsl_odeiv_step_rk4imp; break; } | 1368 | case 6 : {T = gsl_odeiv2_step_rk4imp; break; } |
1369 | case 7 : {T = gsl_odeiv_step_bsimp; break; } | 1369 | case 7 : {T = gsl_odeiv2_step_bsimp; break; } |
1370 | case 8 : {T = gsl_odeiv_step_gear1; break; } | 1370 | case 8 : {T = gsl_odeiv2_step_rk1imp; break; } |
1371 | case 9 : {T = gsl_odeiv_step_gear2; break; } | 1371 | case 9 : {T = gsl_odeiv2_step_msadams; break; } |
1372 | case 10: {T = gsl_odeiv2_step_msbdf; break; } | ||
1372 | default: ERROR(BAD_CODE); | 1373 | default: ERROR(BAD_CODE); |
1373 | } | 1374 | } |
1374 | 1375 | ||
1375 | |||
1376 | gsl_odeiv_step * s = gsl_odeiv_step_alloc (T, xin); | ||
1377 | gsl_odeiv_control * c = gsl_odeiv_control_y_new (eps_abs, eps_rel); | ||
1378 | gsl_odeiv_evolve * e = gsl_odeiv_evolve_alloc (xin); | ||
1379 | |||
1380 | Tode P; | 1376 | Tode P; |
1381 | P.f = f; | 1377 | P.f = f; |
1382 | P.j = jac; | 1378 | P.j = jac; |
1383 | P.n = xin; | 1379 | P.n = xin; |
1384 | 1380 | ||
1385 | gsl_odeiv_system sys = {odefunc, odejac, xin, &P}; | 1381 | gsl_odeiv2_system sys = {odefunc, odejac, xin, &P}; |
1382 | |||
1383 | gsl_odeiv2_driver * d = | ||
1384 | gsl_odeiv2_driver_alloc_y_new (&sys, T, h, eps_abs, eps_rel); | ||
1386 | 1385 | ||
1387 | double t = tsp[0]; | 1386 | double t = tsp[0]; |
1388 | 1387 | ||
1389 | double* y = (double*)calloc(xin,sizeof(double)); | 1388 | double* y = (double*)calloc(xin,sizeof(double)); |
1390 | int i,j; | 1389 | int i,j; |
1390 | int status; | ||
1391 | for(i=0; i< xin; i++) { | 1391 | for(i=0; i< xin; i++) { |
1392 | y[i] = xip[i]; | 1392 | y[i] = xip[i]; |
1393 | solp[i] = xip[i]; | 1393 | solp[i] = xip[i]; |
@@ -1396,22 +1396,24 @@ int ode(int method, double h, double eps_abs, double eps_rel, | |||
1396 | for (i = 1; i < tsn ; i++) | 1396 | for (i = 1; i < tsn ; i++) |
1397 | { | 1397 | { |
1398 | double ti = tsp[i]; | 1398 | double ti = tsp[i]; |
1399 | while (t < ti) | 1399 | |
1400 | { | 1400 | status = gsl_odeiv2_driver_apply (d, &t, ti, y); |
1401 | gsl_odeiv_evolve_apply (e, c, s, | 1401 | |
1402 | &sys, | 1402 | if (status != GSL_SUCCESS) { |
1403 | &t, ti, &h, | 1403 | printf ("error in ode, return value=%d\n", status); |
1404 | y); | 1404 | break; |
1405 | // if (h < hmin) h = hmin; | 1405 | } |
1406 | } | 1406 | |
1407 | // printf ("%.5e %.5e %.5e\n", t, y[0], y[1]); | ||
1408 | |||
1407 | for(j=0; j<xin; j++) { | 1409 | for(j=0; j<xin; j++) { |
1408 | solp[i*xin + j] = y[j]; | 1410 | solp[i*xin + j] = y[j]; |
1409 | } | 1411 | } |
1410 | } | 1412 | } |
1411 | 1413 | ||
1412 | free(y); | 1414 | free(y); |
1413 | gsl_odeiv_evolve_free (e); | 1415 | gsl_odeiv2_driver_free (d); |
1414 | gsl_odeiv_control_free (c); | 1416 | |
1415 | gsl_odeiv_step_free (s); | 1417 | return status; |
1416 | return 0; | ||
1417 | } | 1418 | } |
1419 | |||