diff options
Diffstat (limited to 'lib')
-rw-r--r-- | lib/Numeric/GSL/ODE.hs | 4 | ||||
-rw-r--r-- | lib/Numeric/GSL/gsl-aux.c | 40 |
2 files changed, 22 insertions, 22 deletions
diff --git a/lib/Numeric/GSL/ODE.hs b/lib/Numeric/GSL/ODE.hs index 2251acd..d4f83aa 100644 --- a/lib/Numeric/GSL/ODE.hs +++ b/lib/Numeric/GSL/ODE.hs | |||
@@ -50,8 +50,8 @@ data ODEMethod = RK2 -- ^ Embedded Runge-Kutta (2, 3) method. | |||
50 | | RK2imp -- ^ Implicit 2nd order Runge-Kutta at Gaussian points. | 50 | | RK2imp -- ^ Implicit 2nd order Runge-Kutta at Gaussian points. |
51 | | RK4imp -- ^ Implicit 4th order Runge-Kutta at Gaussian points. | 51 | | RK4imp -- ^ Implicit 4th order Runge-Kutta at Gaussian points. |
52 | | BSimp -- ^ Implicit Bulirsch-Stoer method of Bader and Deuflhard. This algorithm requires the Jacobian. | 52 | | BSimp -- ^ Implicit Bulirsch-Stoer method of Bader and Deuflhard. This algorithm requires the Jacobian. |
53 | | Gear1 -- ^ M=1 implicit Gear method. | 53 | | MSAdams -- ^ A variable-coefficient linear multistep Adams method in Nordsieck form. |
54 | | Gear2 -- ^ M=2 implicit Gear method. | 54 | | MSBDF -- ^ A variable-coefficient linear multistep backward differentiation formula (BDF) method in Nordsieck form. |
55 | deriving (Enum,Eq,Show,Bounded) | 55 | deriving (Enum,Eq,Show,Bounded) |
56 | 56 | ||
57 | -- | A version of 'odeSolveV' with reasonable default parameters and system of equations defined using lists. | 57 | -- | A version of 'odeSolveV' with reasonable default parameters and system of equations defined using lists. |
diff --git a/lib/Numeric/GSL/gsl-aux.c b/lib/Numeric/GSL/gsl-aux.c index 33d7dab..bf3f684 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> |
@@ -1314,33 +1314,33 @@ int ode(int method, double h, double eps_abs, double eps_rel, | |||
1314 | int jac(double, int, const double*, int, int, double*), | 1314 | int jac(double, int, const double*, int, int, double*), |
1315 | KRVEC(xi), KRVEC(ts), RMAT(sol)) { | 1315 | KRVEC(xi), KRVEC(ts), RMAT(sol)) { |
1316 | 1316 | ||
1317 | const gsl_odeiv_step_type * T; | 1317 | const gsl_odeiv2_step_type * T; |
1318 | 1318 | ||
1319 | switch(method) { | 1319 | switch(method) { |
1320 | case 0 : {T = gsl_odeiv_step_rk2; break; } | 1320 | case 0 : {T = gsl_odeiv2_step_rk2; break; } |
1321 | case 1 : {T = gsl_odeiv_step_rk4; break; } | 1321 | case 1 : {T = gsl_odeiv2_step_rk4; break; } |
1322 | case 2 : {T = gsl_odeiv_step_rkf45; break; } | 1322 | case 2 : {T = gsl_odeiv2_step_rkf45; break; } |
1323 | case 3 : {T = gsl_odeiv_step_rkck; break; } | 1323 | case 3 : {T = gsl_odeiv2_step_rkck; break; } |
1324 | case 4 : {T = gsl_odeiv_step_rk8pd; break; } | 1324 | case 4 : {T = gsl_odeiv2_step_rk8pd; break; } |
1325 | case 5 : {T = gsl_odeiv_step_rk2imp; break; } | 1325 | case 5 : {T = gsl_odeiv2_step_rk2imp; break; } |
1326 | case 6 : {T = gsl_odeiv_step_rk4imp; break; } | 1326 | case 6 : {T = gsl_odeiv2_step_rk4imp; break; } |
1327 | case 7 : {T = gsl_odeiv_step_bsimp; break; } | 1327 | case 7 : {T = gsl_odeiv2_step_bsimp; break; } |
1328 | case 8 : {T = gsl_odeiv_step_gear1; break; } | 1328 | case 8 : {T = gsl_odeiv2_step_msadams; break; } |
1329 | case 9 : {T = gsl_odeiv_step_gear2; break; } | 1329 | case 9 : {T = gsl_odeiv2_step_msbdf; break; } |
1330 | default: ERROR(BAD_CODE); | 1330 | default: ERROR(BAD_CODE); |
1331 | } | 1331 | } |
1332 | 1332 | ||
1333 | 1333 | ||
1334 | gsl_odeiv_step * s = gsl_odeiv_step_alloc (T, xin); | 1334 | gsl_odeiv2_step * s = gsl_odeiv2_step_alloc (T, xin); |
1335 | gsl_odeiv_control * c = gsl_odeiv_control_y_new (eps_abs, eps_rel); | 1335 | gsl_odeiv2_control * c = gsl_odeiv2_control_y_new (eps_abs, eps_rel); |
1336 | gsl_odeiv_evolve * e = gsl_odeiv_evolve_alloc (xin); | 1336 | gsl_odeiv2_evolve * e = gsl_odeiv2_evolve_alloc (xin); |
1337 | 1337 | ||
1338 | Tode P; | 1338 | Tode P; |
1339 | P.f = f; | 1339 | P.f = f; |
1340 | P.j = jac; | 1340 | P.j = jac; |
1341 | P.n = xin; | 1341 | P.n = xin; |
1342 | 1342 | ||
1343 | gsl_odeiv_system sys = {odefunc, odejac, xin, &P}; | 1343 | gsl_odeiv2_system sys = {odefunc, odejac, xin, &P}; |
1344 | 1344 | ||
1345 | double t = tsp[0]; | 1345 | double t = tsp[0]; |
1346 | 1346 | ||
@@ -1356,7 +1356,7 @@ int ode(int method, double h, double eps_abs, double eps_rel, | |||
1356 | double ti = tsp[i]; | 1356 | double ti = tsp[i]; |
1357 | while (t < ti) | 1357 | while (t < ti) |
1358 | { | 1358 | { |
1359 | gsl_odeiv_evolve_apply (e, c, s, | 1359 | gsl_odeiv2_evolve_apply (e, c, s, |
1360 | &sys, | 1360 | &sys, |
1361 | &t, ti, &h, | 1361 | &t, ti, &h, |
1362 | y); | 1362 | y); |
@@ -1368,8 +1368,8 @@ int ode(int method, double h, double eps_abs, double eps_rel, | |||
1368 | } | 1368 | } |
1369 | 1369 | ||
1370 | free(y); | 1370 | free(y); |
1371 | gsl_odeiv_evolve_free (e); | 1371 | gsl_odeiv2_evolve_free (e); |
1372 | gsl_odeiv_control_free (c); | 1372 | gsl_odeiv2_control_free (c); |
1373 | gsl_odeiv_step_free (s); | 1373 | gsl_odeiv2_step_free (s); |
1374 | return 0; | 1374 | return 0; |
1375 | } | 1375 | } |