Octave has two built-in functions for solving differential equations. Both are based on reliable ODE solvers written in Fortran.
The function lsode
can be used to solve ODEs of the form
dx -- = f (x, t) dt
using Hindmarsh's ODE solver LSODE.
@anchor{doc-lsode}
dx -- = f(x, t) dt
with
x(t_0) = x_0
The solution is returned in the matrix x, with each row corresponding to an element of the vector t. The first element of t should be @math{t_0} and should correspond to the initial state of the system x_0, so that the first row of the output is x_0.
The first argument, fcn, is a string that names the function to call to compute the vector of right hand sides for the set of equations. The function must have the form
xdot = f (x, t)
in which xdot and x are vectors and t is a scalar.
If fcn is a two-element string array, the first element names the function @math{f} described above, and the second element names a function to compute the Jacobian of @math{f}. The Jacobian function must have the form
jac = j (x, t)
in which jac is the matrix of partial derivatives
| df_1 df_1 df_1 | | ---- ---- ... ---- | | dx_1 dx_2 dx_N | | | | df_2 df_2 df_2 | | ---- ---- ... ---- | df_i | dx_1 dx_2 dx_N | jac = ---- = | | dx_j | . . . . | | . . . . | | . . . . | | | | df_N df_N df_N | | ---- ---- ... ---- | | dx_1 dx_2 dx_N |
The second and third arguments specify the intial state of the system, @math{x_0}, and the initial value of the independent variable @math{t_0}.
The fourth argument is optional, and may be used to specify a set of times that the ODE solver should not integrate past. It is useful for avoiding difficulties with singularities and points where there is a discontinuity in the derivative.
After a successful computation, the value of istate will be 2 (consistent with the Fortran version of LSODE).
If the computation is not successful, istate will be something other than 2 and msg will contain additional information.
You can use the function lsode_options
to set optional
parameters for lsode
.
@anchor{doc-lsode_options}
lsode
.
Given one argument, lsode_options
returns the value of the
corresponding option. If no arguments are supplied, the names of all
the available options and their current values are displayed.
Options include
"absolute tolerance"
"relative tolerance"
abs (local error in x(i)) <= rtol * abs (y(i)) + atol(i)
"integration method"
lsode
will
compute a finite difference approximation of the Jacobian matrix.
"initial step size"
"maximum order"
"maximum step size"
"minimum step size"
"step limit"
Here is an example of solving a set of three differential equations using
lsode
. Given the function
function xdot = f (x, t) xdot = zeros (3,1); xdot(1) = 77.27 * (x(2) - x(1)*x(2) + x(1) \ - 8.375e-06*x(1)^2); xdot(2) = (x(3) - x(1)*x(2) - x(2)) / 77.27; xdot(3) = 0.161*(x(1) - x(3)); endfunction
and the initial condition x0 = [ 4; 1.1; 4 ]
, the set of
equations can be integrated using the command
t = linspace (0, 500, 1000); y = lsode ("f", x0, t);
If you try this, you will see that the value of the result changes dramatically between t = 0 and 5, and again around t = 305. A more efficient set of output points might be
t = [0, logspace (-1, log10(303), 150), \ logspace (log10(304), log10(500), 150)];
Octave also includes ODESSA, a modification of LSODE that allows for the computation of parametric sensitivity information simultaneously with the solution of the set of ODEs.
@anchor{doc-odessa}
dx -- = f(x, t; p) dt
with
x(t_0) = x_0
and simultaneously compute the first-order sensitivity coefficients given by
s'(t) = j(t)*s(t) + df/dp
in which
s(t) = dx(t)/dp (sensitivity functions) s'(t) = d(dx(t)/dp)/dt j(t) = df(x,t;p)/dx(t) (Jacobian matrix) df/dp = df(x,t;p)/dp (inhomogeneity matrix)
The solution is returned in the matrix x, with each row corresponding to an element of the vector t. The first element of t should be @math{t_0} and should correspond to the initial state of the system x_0, so that the first row of the output is x_0.
The sensitivities are returned in a list of matrices, sx, with each element of the list corresponding to an element of the vector t.
The first argument, fcn, is a string that names the function to call to compute the vector of right hand sides for the set of equations. The function must have the form
xdot = f (x, t, p)
in which xdot and x are vectors and t is a scalar.
The variable p is a vector of parameters.
The fcn argument may also be an array of strings
["f"; "j"; "b"]
in which the first element names the function @math{f} described above, the second element names a function to compute the Jacobian of @math{f}, and the third element names a function to compute the inhomogeneity matrix.
The Jacobian function must have the form
jac = j (x, t, p)
in which x, t, and p have the same meanings as above for the function f, and jac is the matrix of partial derivatives
df_i jac = ---- dx_j
The function b must have the form
dfdp = b (x, t, p, c)
in which x, t, and p have the same meanings as above for the function f, c indicates which partial derivatives to return in dfdp. For example, if c is 2, you should return the vector
df_i dfdp = ----, i = 1:length(x) dp_2
The second argument, x_0, specifies the intial state of the system.
The third argument, p, specifies the set of parameters.
The fourth argument, sx_0 specifies the initial values of the sensitivities.
The sixth argument is optional, and may be used to specify a set of times that the ODE solver should not integrate past. It is useful for avoiding difficulties with singularities and points where there is a discontinuity in the derivative.
After a successful computation, the value of istate will be 2 (consistent with the Fortran version of ODESSA).
If the computation is not successful, istate will be something other than 2 and msg will contain additional information.
You can use the function lsode_options
to set optional
parameters for lsode
.
@anchor{doc-odessa_options}
odessa
.
Given one argument, odessa_options
returns the value of the
corresponding option. If no arguments are supplied, the names of all
the available options and their current values are displayed.
Options include
"absolute tolerance"
"relative tolerance"
abs (local error in x(i)) <= rtol * abs (y(i)) + atol(i)
"integration method"
lsode
will
compute a finite difference approximation of the Jacobian matrix.
"initial step size"
"maximum order"
"maximum step size"
"minimum step size"
"step limit"
See Alan C. Hindmarsh, ODEPACK, A Systematized Collection of ODE
Solvers, in Scientific Computing, R. S. Stepleman, editor, (1983) for
more information about the inner workings of lsode
.
The function daspk
can be used to solve DAEs of the form
0 = f (x-dot, x, t), x(t=0) = x_0, x-dot(t=0) = x-dot_0
using Petzold's DAE solver DASPK.
@anchor{doc-daspk}
0 = f (xdot, x, t)
with
x(t_0) = x_0, xdot(t_0) = xdot_0
The solution is returned in the matrices x and xdot, with each row in the result matrices corresponding to one of the elements in the vector t. The first element of t should be @math{t_0} and correspond to the initial state of the system x_0 and its derivative xdot_0, so that the first row of the output x is x_0 and the first row of the output xdot is xdot_0.
The first argument, fcn, is a string that names the function to call to compute the vector of residuals for the set of equations. It must have the form
res = f (x, xdot, t)
in which x, xdot, and res are vectors, and t is a scalar.
If fcn is a two-element string array, the first element names the function @math{f} described above, and the second element names a function to compute the modified Jacobian
df df jac = -- + c ------ dx d xdot
The modified Jacobian function must have the form
jac = j (x, xdot, t, c)
The second and third arguments to daspk
specify the initial
condition of the states and their derivatives, and the fourth argument
specifies a vector of output times at which the solution is desired,
including the time corresponding to the initial condition.
The set of initial states and derivatives are not strictly required to
be consistent. If they are not consistent, you must use the
daspk_options
function to provide additional information so
that daspk
can compute a consistent starting point.
The fifth argument is optional, and may be used to specify a set of times that the DAE solver should not integrate past. It is useful for avoiding difficulties with singularities and points where there is a discontinuity in the derivative.
After a successful computation, the value of istate will be greater than zero (consistent with the Fortran version of DASPK).
If the computation is not successful, the value of istate will be less than zero and msg will contain additional information.
You can use the function daspk_options
to set optional
parameters for daspk
.
@anchor{doc-daspk_options}
daspk
.
Given one argument, daspk_options
returns the value of the
corresponding option. If no arguments are supplied, the names of all
the available options and their current values are displayed.
Options include
"absolute tolerance"
"relative tolerance"
abs (local error in x(i)) <= rtol(i) * abs (Y(i)) + atol(i)
"compute consistent initial condition"
ddaspk
can solve
one of two initialization problems:
"algebraic variables"
option to declare which variables in the
problem are algebraic.
"use initial condition heuristics"
"initial condition heuristics"
MXNIT
MXNJ
MXNH
"compute consistent initial condition"
option has
been set to 1 (default is 5).
Note that the maximum number of Newton iterations allowed in all is
MXNIT*MXNJ*MXNH
if the "compute consistent initial
condition"
option has been set to 1 and MXNIT*MXNJ
if it is
set to 2.
LSOFF
STPTOL
EPINIT
EPINIT*EPCON
,
where EPCON
= 0.33 is the analogous test constant used in the
time steps. The default is EPINIT
= 0.01.
"print initial condition info"
"exclude algebraic variables from error test"
"algebraic variables"
option to
declare which variables in the problem are algebraic (default is 0).
"algebraic variables"
compute consistent initial condition"
and
"exclude algebraic variables from error test"
options.
"enforce inequality constraints"
"inequality constraint types"
option (default is 0).
"inequality constraint types"
"enforce inequality constraints"
option is nonzero.
"initial step size"
"maximum order"
"maximum step size"
Octave also includes DASSL, an earlier version of Daspk, and dasrt, which can be used to solve DAEs with constraints (stopping conditions).
@anchor{doc-dasrt}
0 = f (xdot, x, t)
with
x(t_0) = x_0, xdot(t_0) = xdot_0
with functional stopping criteria (root solving).
The solution is returned in the matrices x and xdot, with each row in the result matrices corresponding to one of the elements in the vector t_out. The first element of t should be @math{t_0} and correspond to the initial state of the system x_0 and its derivative xdot_0, so that the first row of the output x is x_0 and the first row of the output xdot is xdot_0.
The vector t provides an upper limit on the length of the integration. If the stopping condition is met, the vector t_out will be shorter than t, and the final element of t_out will be the point at which the stopping condition was met, and may not correspond to any element of the vector t.
The first argument, fcn, is a string that names the function to call to compute the vector of residuals for the set of equations. It must have the form
res = f (x, xdot, t)
in which x, xdot, and res are vectors, and t is a scalar.
If fcn is a two-element string array, the first element names the function @math{f} described above, and the second element names a function to compute the modified Jacobian
df df jac = -- + c ------ dx d xdot
The modified Jacobian function must have the form
jac = j (x, xdot, t, c)
The optional second argument names a function that defines the constraint functions whose roots are desired during the integration. This function must have the form
g_out = g (x, t)
and return a vector of the constraint function values. If the value of any of the constraint functions changes sign, DASRT will attempt to stop the integration at the point of the sign change.
If the name of the constraint function is omitted, dasrt
solves
the saem problem as daspk
or dassl
.
Note that because of numerical errors in the constraint functions due to roundoff and integration error, DASRT may return false roots, or return the same root at two or more nearly equal values of T. If such false roots are suspected, the user should consider smaller error tolerances or higher precision in the evaluation of the constraint functions.
If a root of some constraint function defines the end of the problem, the input to DASRT should nevertheless allow integration to a point slightly past that root, so that DASRT can locate the root by interpolation.
The third and fourth arguments to dasrt
specify the initial
condition of the states and their derivatives, and the fourth argument
specifies a vector of output times at which the solution is desired,
including the time corresponding to the initial condition.
The set of initial states and derivatives are not strictly required to be consistent. In practice, however, DASSL is not very good at determining a consistent set for you, so it is best if you ensure that the initial values result in the function evaluating to zero.
The sixth argument is optional, and may be used to specify a set of times that the DAE solver should not integrate past. It is useful for avoiding difficulties with singularities and points where there is a discontinuity in the derivative.
After a successful computation, the value of istate will be greater than zero (consistent with the Fortran version of DASSL).
If the computation is not successful, the value of istate will be less than zero and msg will contain additional information.
You can use the function dasrt_options
to set optional
parameters for dasrt
.
@anchor{doc-dasrt_options}
dasrt
.
Given one argument, dasrt_options
returns the value of the
corresponding option. If no arguments are supplied, the names of all
the available options and their current values are displayed.
Options include
"absolute tolerance"
"relative tolerance"
abs (local error in x(i)) <= rtol(i) * abs (Y(i)) + atol(i)
"initial step size"
"maximum order"
"maximum step size"
"step limit"
See K. E. Brenan, et al., Numerical Solution of Initial-Value Problems in Differential-Algebraic Equations, North-Holland (1989) for more information about the implementation of DASSL.
Go to the first, previous, next, last section, table of contents.