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_odeiv_step_name (s));would print something like
step method is 'rk4'
.
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 applies the stepping function s to the system of equations defined by dydt, 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.
If the user-supplied functions defined in the system dydt return a status other than
GSL_SUCCESS
the step will be aborted. In this 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. To distinguish between error codes from the user-supplied functions and those fromgsl_odeiv_step_apply
itself, any user-defined return values should be distinct from the standard GSL error codes.
The following algorithms are available,
Embedded Runge-Kutta-Fehlberg (4, 5) method. This method is a good general-purpose integrator.