view scripts/ode/odeplot.m @ 31263:449ed6f427cb

ode45/23/23s: Implement Events, OutputFcn, & Refine options (bug #49408 and #63063) * scripts/ode/ode23.m: Remove disabling of Refine option with struct output. Modify solution struct to output two sets of solution variables: output_t, output_x and ode_t and ode_x, and transpose struct output variables for improved Matlab compatibility. Update BISTs and perform minor code formatting. * scripts/ode/ode23s.m: Make same changes as ode23.m. * scripts/ode/ode45.m: Make same changes as ode23.m. Remove comment indicating that Refine is not implemented. * scripts/ode/private/integrate_adaptive.m: Update internal handling of variables t and x, separating them into ode_t & ode_x for internal integration and output_t & output_x for function output or calls to OutputFcn. Replace prior attempt at Refine option with new implementation. Specify time output or Refine != 0 are both interpolated from internal variables (ode_t, ode_x) for output of non-struct variables and/or for use with OutputFcn. Improve event handling when multiple Events (including at least one terminal Event) are detected in a single simulation step so that all Events up to and including the first terminal one are reported, and final data point is set to that of terminal Event. Send multiple data points in a single call to OutputFcn if they are all interpolated from a single integration step. Remove warning for termination when term signal is received from Events or OutputFcn. Return both internal variables (ode_t, ode_x) and interpolated variables (output_t, output_x) to allow calling function to correctly return either struct or separate variables. * scripts/ode/private/ode_event_handler.m: Sort multiple Events in ascending time order when they are encountered in one integration step. Remove any events after the time of a terminal Event. * scripts/ode/odeset.m: Update docstring to remove indication that Refine is not implemented * scripts/ode/odeplot.m: Update docstring to indicate that input t can be a scalar or vector. Add file test. * etc/NEWS.8.md: Add descriptions of changes under General improvements and Matlab compatibility.
author Ken Marek <marek_ka@mercer.edu>
date Wed, 05 Oct 2022 16:53:01 -0400
parents 796f54d4ddbf
children c332a2f2959f
line wrap: on
line source

########################################################################
##
## Copyright (C) 2006-2022 The Octave Project Developers
##
## See the file COPYRIGHT.md in the top-level directory of this
## distribution or <https://octave.org/copyright/>.
##
## This file is part of Octave.
##
## Octave is free software: you can redistribute it and/or modify it
## under the terms of the GNU General Public License as published by
## the Free Software Foundation, either version 3 of the License, or
## (at your option) any later version.
##
## Octave is distributed in the hope that it will be useful, but
## WITHOUT ANY WARRANTY; without even the implied warranty of
## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
## GNU General Public License for more details.
##
## You should have received a copy of the GNU General Public License
## along with Octave; see the file COPYING.  If not, see
## <https://www.gnu.org/licenses/>.
##
########################################################################

## -*- texinfo -*-
## @deftypefn {} {@var{stop_solve} =} odeplot (@var{t}, @var{y}, @var{flag})
##
## Open a new figure window and plot the solution of an ode problem at each
## time step during the integration.
##
## The types and values of the input parameters @var{t} and @var{y} depend on
## the input @var{flag} that is of type string.  Valid values of @var{flag}
## are:
##
## @table @option
## @item @qcode{"init"}
## The input @var{t} must be a column vector of length 2 with the first and
## last time step (@code{[@var{tfirst} @var{tlast}]}.  The input @var{y}
## contains the initial conditions for the ode problem (@var{y0}).
##
## @item @qcode{""}
## The input @var{t} must be a scalar double or vector specifying the time(s)
## for which the solution in input @var{y} was calculated.
##
## @item @qcode{"done"}
## The inputs should be empty, but are ignored if they are present.
## @end table
##
## @code{odeplot} always returns false, i.e., don't stop the ode solver.
##
## Example: solve an anonymous implementation of the
## @nospell{@qcode{"Van der Pol"}} equation and display the results while
## solving.
##
## @example
## @group
## fvdp = @@(t,y) [y(2); (1 - y(1)^2) * y(2) - y(1)];
##
## opt = odeset ("OutputFcn", @@odeplot, "RelTol", 1e-6);
## sol = ode45 (fvdp, [0 20], [2 0], opt);
## @end group
## @end example
##
## Background Information:
## This function is called by an ode solver function if it was specified in
## the @qcode{"OutputFcn"} property of an options structure created with
## @code{odeset}.  The ode solver will initially call the function with the
## syntax @code{odeplot ([@var{tfirst}, @var{tlast}], @var{y0}, "init")}.  The
## function initializes internal variables, creates a new figure window, and
## sets the x limits of the plot.  Subsequently, at each time step during the
## integration the ode solver calls @code{odeplot (@var{t}, @var{y}, [])}.
## At the end of the solution the ode solver calls
## @code{odeplot ([], [], "done")} so that odeplot can perform any clean-up
## actions required.
## @seealso{odeset, odeget, ode23, ode45}
## @end deftypefn

function stop_solve = odeplot (t, y, flag)

  ## No input argument checking is done for better performance
  persistent hlines num_lines told yold;

  ## odeplot never stops the integration
  stop_solve = false;

  if (isempty (flag))
    ## Default case, plot and return a value
    told = [told; t(:)];
    yold = [yold, y];
    for i = 1:num_lines
      set (hlines(i), "xdata", told, "ydata", yold(i,:));
    endfor
    drawnow ();

    retval = false;

  elseif (strcmp (flag, "init"))
    ## t is either the time slot [tstart tstop] or [t0, t1, ..., tn]
    ## y is the initial value vector for the ode solution
    told = t(1);
    yold = y(:);
    figure ();
    hlines = plot (told, yold, "o-");
    xlim ([t(1), t(end)]);  # Fix limits which also speeds up plotting
    num_lines = numel (hlines);

  elseif (strcmp (flag, "done"))
    ## Cleanup after ode solver has finished.
    hlines = num_lines = told = yold = [];

  endif

endfunction


%!demo
%! ## Solve an anonymous implementation of the Van der Pol equation
%! ## and display the results while solving
%! fvdp = @(t,y) [y(2); (1 - y(1)^2) * y(2) - y(1)];
%! opt = odeset ("OutputFcn", @odeplot, "RelTol", 1e-6);
%! sol = ode45 (fvdp, [0 20], [2 0], opt);

## FIXME: convert to demo or a visible=off test with failable assert/error
##        statemments
## Test that the function works for expected ode OutputFcn calls
## %!test
## %! t = linspace(0,2,10);
## %! y = [0.2*t; -0.1*t.^2-1; sin(t)];
## %! odeplot ([0 2], y(:,1), "init");
## %! odeplot (t(2), y(:,2), []);
## %! odeplot (t(3:end), y(:, 3:end), []);
## %! odeplot ([], [], "done");
## %! close all