octave: Ordinary Differential Equations

 
 24.1 Ordinary Differential Equations
 ====================================
 
 The function ‘lsode’ can be used to solve ODEs of the form
 
      dx
      -- = f (x, t)
      dt
 
 using Hindmarsh’s ODE solver LSODE.
 
  -- : [X, ISTATE, MSG] = lsode (FCN, X_0, T)
  -- : [X, ISTATE, MSG] = lsode (FCN, X_0, T, T_CRIT)
      Ordinary Differential Equation (ODE) solver.
 
      The set of differential equations to solve is
 
           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 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, inline, or function handle
      that names the function f 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 or a two-element cell array of
      strings, inline functions, or function handles, the first element
      names the function f described above, and the second element names
      a function to compute the Jacobian of 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 argument specifies the initial state of the system x_0.
      The third argument is a vector, T, specifying the time values for
      which a solution is sought.
 
      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’.
 
DONTPRINTYET       See also: Seedaspk XREFdaspk, Seedassl XREFdassl, *noteDONTPRINTYET       See also: Seedaspk XREFdaspk, Seedassl XREFdassl, See
      dasrt XREFdasrt.
 
  -- : lsode_options ()
  -- : val = lsode_options (OPT)
  -- : lsode_options (OPT, VAL)
      Query or set options for the function ‘lsode’.
 
      When called with no arguments, the names of all available options
      and their current values are displayed.
 
      Given one argument, return the value of the option OPT.
 
      When called with two arguments, ‘lsode_options’ sets the option OPT
      to value VAL.
 
      Options include
 
      "absolute tolerance"
           Absolute tolerance.  May be either vector or scalar.  If a
           vector, it must match the dimension of the state vector.
 
      "relative tolerance"
           Relative tolerance parameter.  Unlike the absolute tolerance,
           this parameter may only be a scalar.
 
           The local error test applied at each integration step is
 
                  abs (local error in x(i)) <= ...
                      rtol * abs (y(i)) + atol(i)
 
      "integration method"
           A string specifying the method of integration to use to solve
           the ODE system.  Valid values are
 
           "adams"
           "non-stiff"
                No Jacobian used (even if it is available).
 
           "bdf"
           "stiff"
                Use stiff backward differentiation formula (BDF) method.
                If a function to compute the Jacobian is not supplied,
                ‘lsode’ will compute a finite difference approximation of
                the Jacobian matrix.
 
      "initial step size"
           The step size to be attempted on the first step (default is
           determined automatically).
 
      "maximum order"
           Restrict the maximum order of the solution method.  If using
           the Adams method, this option must be between 1 and 12.
           Otherwise, it must be between 1 and 5, inclusive.
 
      "maximum step size"
           Setting the maximum stepsize will avoid passing over very
           large regions (default is not specified).
 
      "minimum step size"
           The minimum absolute step size allowed (default is 0).
 
      "step limit"
           Maximum number of steps allowed (default is 100000).
 
    Here is an example of solving a set of three differential equations
 using ‘lsode’.  Given the function
 
      ## oregonator differential equation
      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)];
 
    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’.
 
    An m-file for the differential equation used above is included with
 the Octave distribution in the examples directory under the name
 ‘oregonator.m’.
 

Menu