Sundials ordinary differential equation solver
[t,y] = cvode(f,tspan,y0,options) [t,y,te,ye,ie] = cvode(f,tspan,y0,options) sol = cvode(...) solext = cvode(sol,tfinal,options)
f | a function, a string or a list, the right hand side of the differential equation. |
tspan | double vector, time interval or time points. |
y0 | a double array: initial state of the ode. |
options | a sequence of optional named arguments (see sundials options) |
t | vector of time points used by the solver. |
y | array of solution at time values in t |
te,ye,ie | time of event, solution and index of event. |
sol, solext | MList of |
tfinal | final time of extended solution |
cvode
computes the solution of real or complex ordinary different equations defined by . It is an interface to CVODE solver of Sundials library with default Adams method (BDF method can be chosen in options). The simplest call of cvode is when
tspan
is a two component vector:
[t,y] = cvode(f,[t0 tf],y0)
where y0
is the
initial condition, t0, tf
are the initial and final time, y
is the
array of solutions [y(t(1)),y(t(2)),...]
at time values in t, where concatenation is done on next dimension if y0
is not a vector. The time values in t
are those used by the solver to meet default relative and absolute estimated local error (which can be changed in options). In the simplest case the right hand side is computed by a Scilab function with two arguments, for example for the ode the right hand side is coded as
function yprime = f(t,y) yprime = t+y end
See the CVODE user functions to learn how to pass extra parameters and/or use DLL entrypoints. When t
has more than two components:
[t,y] = cvode(f,[t0 t1 ... tf],y0)
the solution is computed at prescribed points with the same precision as the two components syntax. However, the result may slightly differ (within chosen tolerance) since t1-t0
may give a different guess of the initial step used by the solver. Hence, solver internal steps are roughly the same and solution at user points is computed by continuous extension formulas.
When searching for where some functions defined in
options
vanish (see the Events section below) the syntax
[t,y,te,ye,ie] = cvode(f,tspan,y0,options)
allows to recover values in
te,ye
and in ie
the number(s) of the vanishing function.
When solution has to be further evaluated at arbitrary points which are not known in advance, then the syntax
sol = cvode(f,[t0 tf],y0)
yields an MList of _odeSolution
type, which can be used later as an interpolant, see the solution output section.
In the following example, we solve the Ordinary Differential Equation
with the initial
condition
and
. Since the original equation is of second order
we define the state of the equivalent first order equation as
. The right hand side is computed by the function
%cvode_vdp1
(defined in the Sundials module):
function vdot=%cvode_vdp1(t, v) mu = 1; vdot = [v(2); mu*(1-v(1)^2)*v(2)-v(1)] end
[t,v] = cvode(%cvode_vdp1, [0 10], [0; 2]); plot(t,v) | ![]() | ![]() |
Agebraic equations to be solved simultaneously with the integration of the ODE can be specified by using the events
option:
[t,y,te,ye,ie] = cvode(f,tspan,y0,events = g)
where g
in its simplest form is a Scilab function with prototype [eq,term,dir] = g(t,y)
, where eq(i)=0
when event i
occurs, term
is a boolean vector, with term(i)=%t
value if integration has to be stopped when event i
occurs. Vector dir
allows to
select event direction, with term(i)
can take the values -1,1
if solution has to be decreasing or increasing or 0
if direction does not matter. If the corresponding
behavior does not matter, dir
or both term,dir
outputs can be omitted in the function prototype. For example, in order to find all points in [0,10] where the solution of Van Der Pol equation veryfies
, we use the following code:
function eq=g(t, v) eq = v(1)-1; end [t,v,te,ve,ie] = cvode(%cvode_vdp1, [0 10], [0; 2], events = g); plot(te,ve(1,:),"or",t,v) | ![]() | ![]() |
The syntax sol = cvode(f,[t0 tf],y0)
yields an MList with fields sol.solver
(the solver name, here "cvode"), sol.method
("adams" or "bdf"), sol.t
(the solver time steps), sol.y
(the solution at solver timesteps), sol.idata
(a pointer to an internal CVODE object) and sol.stats
, a structure gathering the solver statistics. When events have been considered additional fields sol.te, sol.ye, sol.ie
are also present.
For an arbitrary time vector t
, with values in the interval [t0,tf]
, y=sol(t)
yields by costless interpolation the same approximation as if components of t
where used in tspan
. The MList sol
can also be used to restart the solver and extend the solution for t>tf
(see the following section). The derivative of solution can also be obtained with [y,yd]=sol(t)
and a particular component i
of the solution is obtained with y=sol(t,i)
. For example, the solution of the Van Der Pol equation can be refined in [3,5] with the following code
The syntax solext = cvode(sol,tfinal)
extends the solution by restarting the solver from initial time sol.t($)
and initial condition sol.y(:,$)
. It can be used when you know in advance that the solution is not differentiable at some time point. By stopping then restarting the solver at time of discontinuity you optimize solver effort by avoiding very small time steps. In the options
sequence you can change almost all previously used options used in the call of cvode
which yielded sol
. The right hand side and the initial condition can be overriden by using the specific f
option and y0
option, respectively.
sol = cvode(%cvode_vdp1, [0 5], [0; 2]); solext = cvode(sol, 10); t = linspace(0,5,500); plot(t,sol(t,1)) t = linspace(5,10,500); plot(t,solext(t,1),'r') | ![]() | ![]() |
// Complex ode function out=crhs(t, y) out = 10*exp(2*%i*%pi*t)*y; endfunction [t,z] = cvode(crhs,[0 5],1); plot(t,real(z),t,imag(z)) legend "real(z)" "imag(z)" | ![]() | ![]() |
A. C. Hindmarsh, P. N. Brown, K. E. Grant, S. L. Lee, R. Serban, D. E. Shumaker, and C. S. Woodward, "SUNDIALS: Suite of Nonlinear and Differential/Algebraic Equation Solvers," ACM Transactions on Mathematical Software, 31(3), pp. 363-396, 2005. Also available as LLNL technical report UCRL-JP-200037.