Interactive trajectory evolution

Without timeseries

InteractiveDynamics.interactive_evolutionFunction
interactive_evolution(ds::DynamicalSystem, u0s; kwargs...)

Launch an interactive application that can evolve the initial conditions u0s (vector of vectors) of the given dynamical system. All initial conditions are evolved in parallel and at exactly the same time. Two controls allow you to pause/resume the evolution and to adjust the speed. The application can run forever (trajectories are computed on demand).

The function returns fig, obs. fig is the overarching fig (the entire GUI) and can be recorded. obs is a vector of observables, each containing the current state of the trajectory.

Keywords

  • transform = identity : Transformation applied to the state of the dynamical system before plotting. Can even return a vector that is of higher dimension than ds.
  • idxs = 1:min(length(transform(ds.u0)), 3) : Which variables to plot (up to three can be chosen). Variables are selected after transform has been applied.
  • colors : The color for each trajectory. Random colors are chosen by default.
  • lims : A tuple of tuples (min, max) for the axis limits. If not given, they are automatically deduced by evolving each of u0s 1000 units and picking most extreme values (limits cannot be adjust after application is launched).
  • m = 1.0 : The trajectory endpoints have a marker. A heuristic is done to choose appropriate marker size given the trajectory size. m is a multiplier that scales the marker size.
  • plotkwargs = NamedTuple() : A named tuple of keyword arguments propagated to the plotting function (lines for continuous, scatter for discrete systems). plotkwargs can also be a vector of named tuples, in which case each initial condition gets different arguments.
  • tail = 1000 : Length of plotted trajectory (in step units).
  • diffeq = DynamicalSystems.CDS_KWARGS : Named tuple of keyword arguments propagated to the solvers of DifferentialEquations.jl (for continuous systems). Because trajectories are not pre-computed and interpolated, it is recommended to use a combination of arguments that limit maximum stepsize, to ensure smooth curves. For example:
    using OrdinaryDiffEq
    diffeq = (alg = Tsit5(), dtmax = 0.01)
source

To generate the video at the start of this page run

using InteractiveDynamics
using DynamicalSystems, GLMakie
using OrdinaryDiffEq

ds = Systems.henonheiles()  # 4D chaotic/regular continuous system

u0s = [[0.0, -0.25, 0.42081, 0.0],
[0.0, 0.1, 0.5, 0.0],
[0.0, -0.31596, 0.354461, 0.0591255]]

diffeq = (alg = Vern9(), dtmax = 0.01)
idxs = (1, 2, 4)
colors = ["#233B43", "#499cbf", "#E84646"]

figure, obs = interactive_evolution(
    ds, u0s; idxs, tail = 10000, diffeq, colors
)

And here is another version for a discrete system:

using InteractiveDynamics
using DynamicalSystems, GLMakie

ds = Systems.towel() # 3D chaotic discrete system
u0s = [0.1ones(3) .+ 1e-3rand(3) for _ in 1:3]

figure, obs = interactive_evolution(
    ds, u0s; idxs = SVector(1, 2, 3), tail = 100000,
)

With timeseries

InteractiveDynamics.interactive_evolution_timeseriesFunction
interactive_evolution_timeseries(ds::DynamicalSystem, u0s, ps = nothing; kwargs...)

If ps === nothing, this function does the same as interactive_evolution, but in addition to the state space plot there is a panel with the variable timeseries plotted and animated in real time.

If ps is not nothing, then it must be a dictionary, mapping keys of the system parameter container (ds.p) to possible ranges of values. The app then will add some additional controls on the left side which allow one to interactively change system parameters and then click the "update" button to translate the new parameters to system evolution. This can be done without stopping the live system evolution. Notice that in this scenario it is recommended to provide the lims keyword manually.

The following additional keywords also apply:

  • total_span : How much the x-axis of the timeseries plots should span (in real time units)
  • linekwargs = NamedTuple() : Extra keywords propagated to the timeseries plots.
  • pnames = Dict(keys(ps) .=> keys(ps)) : Dictionary mapping parameter keys to labels. Only valid if ps is a dictionary and not nothing.
source