Next: Adaptive Step-size Control, Previous: Defining the ODE System, Up: Ordinary Differential Equations
The lowest level components are the stepping functions which advance a solution from time t to t+h for a fixed step-size h and estimate the resulting local error.
This function returns a pointer to a newly allocated instance of a stepping function of type T for a system of dim dimensions.
This function resets the stepping function s. It should be used whenever the next use of s will not be a continuation of a previous step.
This function frees all the memory associated with the stepping function s.
This function returns a pointer to the name of the stepping function. For example,
printf ("step method is '%s'\n", gsl_odeiv2_step_name (s));would print something like
step method is 'rkf45'
.
This function returns the order of the stepping function on the previous step. This order can vary if the stepping function itself is adaptive.
This function can be called after allocation of stepper s and control object c to allow the stepper to obtain desired error level via the control object. This is a requirement for some steppers.
This function applies the stepping function s to the system of equations defined by sys, using the step-size h to advance the system from time t and state y to time t+h. The new state of the system is stored in y on output, with an estimate of the absolute error in each component stored in yerr. If the argument dydt_in is not null it should point an array containing the derivatives for the system at time t on input. This is optional as the derivatives will be computed internally if they are not provided, but allows the reuse of existing derivative information. On output the new derivatives of the system at time t+h will be stored in dydt_out if it is not null.
The stepping function returns
GSL_FAILURE
if it is unable to compute the requested step. Also, if the user-supplied functions defined in the system sys return a status other thanGSL_SUCCESS
the step will be aborted. In that case, the elements of y will be restored to their pre-step values and the error code from the user-supplied function will be returned. Failure may be due to a singularity in the system or too large step-size h. In that case the step should be attempted again with a smaller step-size, e.g. h/2.If the control object is not appropriately set via
gsl_odeiv2_step_set_control
for those steppers that need it, the stepping function returnsGSL_EFAULT
.
The following algorithms are available,
Explicit 4th order (classical) Runge-Kutta. Error estimation is carried out by the step doubling method. For more efficient estimate of the error, use the explicit methods described below.
Explicit embedded Runge-Kutta-Fehlberg (4, 5) method. This method is a good general-purpose integrator.
Implicit Gaussian first order Runge-Kutta. Also known as implicit Euler or backward Euler method. Error estimation is carried out by the step doubling method. This algorithm requires the Jacobian and the setting of the control object via
gsl_odeiv2_step_set_control
.
Implicit Gaussian second order Runge-Kutta. Also known as implicit mid-point rule. Error estimation is carried out by the step doubling method. This stepper requires the Jacobian and the setting of the control object via
gsl_odeiv2_step_set_control
.
Implicit Gaussian 4th order Runge-Kutta. Error estimation is carried out by the step doubling method. This algorithm requires the Jacobian and the setting of the control object via
gsl_odeiv2_step_set_control
.
Implicit Bulirsch-Stoer method of Bader and Deuflhard. The method is generally suitable for stiff problems. This stepper requires the Jacobian.
A variable-coefficient linear multistep Adams method in Nordsieck form. This stepper uses explicit Adams-Bashforth (predictor) and implicit Adams-Moulton (corrector) methods in P(EC)^m functional iteration mode. Method order varies dynamically between 1 and 12. This stepper requires the setting of the control object via
gsl_odeiv2_step_set_control
.
A variable-coefficient linear multistep backward differentiation formula (BDF) method in Nordsieck form. This stepper uses the explicit BDF formula as predictor and implicit BDF formula as corrector. A modified Newton iteration method is used to solve the system of non-linear equations. Method order varies dynamically between 1 and 5. The method is generally suitable for stiff problems. This stepper requires the Jacobian and the setting of the control object via
gsl_odeiv2_step_set_control
.