API
This page is a dump of all the docstrings found in the code.
ODEHybrid.ODE_SET
— TypeCopy of MATLAB's ODE SET struct
ODEHybrid.TimeSeriesLogger
— TypeA struct for keeping track of numerous pieces of data over time throughout a process. It's designed to be easy to use and relatively quickly. It logs numeric types of any size (scalars, vectors, or matrices), as long as the size is consistent from sample to sample.
Example:
log = TimeSeriesLogger(); # Make a new logger.
x = [0; 0];
for t = 1:100 # Make a random walk.
x = x + randn(2, 1); # Take a single step.
add!(log, "walk", t, x); # Log a single sample, x, at time t.
end
tsl_plot(log) # Plot everything.
# We can also access specific logs by their names. In the above, we only
# have one log ('walk'). Let's get that log from the logger.
t_log, x_log = get_log(log, "walk"); # Get a specific log.
plot(x_log[:, 1], x_log[:, 2]); # Do something with it.
# We can make sure a log exists before trying to do anything with it:
if contains(log,"walk")
x_log = get_log(log, "walk");
plot(x_log[:, 1], x_log[:, 2]);
;
# If we want to log something but don't want it plotted when we call
# |plot|, then we can pass in false to the logger when we add the signal.
add!(log, "var1", t, x, false);
# To group items into different figures when they're plotted, give them
# "groups":
add!(log, "var1", t, x, true, "group1");
add!(log, "var2", t, y, true, "group2");
add!(log, "foo", t, bar, true, "group1");
# All of the items with common group names will be added to the same
# figures as subplots.
# Finally, we can clear out the logger with |initialize|. This deletes all
# data and returns the logger to its initial state.
initialize!(log);
Methods
- initialize! Clear out all data.
- add! Add a single data point to a time series.
- contains See if a time series exists (by name).
- get_log Return a log by name.
- tsl_plot Plot the logged series.
Notes
TimeSeriesLogger was created as the logging mechanism for odehybrid. However, it's not dependent on that project and so can be used anywhere.
TimeSeriesLogger is not related to the 'timeseries' class in MATLAB.
Since this class never knows how much more data is going to be logged, it can't preallocate the appropriate amount of space. However, increasing the size of its data stores by one row every time a new sample is added is very slow. To combat this, the data store starts off small and doubles whenever the store is at capacity. Say we're logging a 4x1 signal. When the first sample is added (say it's log k), data{k}{2} will be 1x4. When the second signal is added, it becomes 2x4. For the third, it's 4x4, then 8x4, 16x4, etc. A separate counter stores how much of this allocated space is currently used. This reduces the number of allocations from n to log2(n). Practically, it saves a little time during logging without too much complexity.
Online doc: http://www.anuncommonlab.com/doc/odehybrid/TimeSeriesLogger.html
Copyright 2014 An Uncommon Lab
ODEHybrid.add!
— FunctionAdd a new log or append to an existing log by name. Returns true iff a new log was created.
add!(log, "var1", t, x);
add!(log, "var1", t, x, true); # Show log in plots (default)
add!(log, "var1", t, x, false); # Don't show log in plots
# Signals can be grouped together into figures by given them
# a common group argument. Here, both var1 and var2 logs will
# be plotted together.
add!(log, "var1", t, x, true, "group1");
add!(log, "var2", t2, y, true, "group1");
NOTE that unlike in MATLAB, the signals are stored separately in times (t) and data (x), (e.g., the kth pair would be log.times[k] and log.data[k])
ODEHybrid.contains
— MethodReturn 'true' iff the TimeSeriesLogger contains this name
ODEHybrid.get_log
— MethodReturn a specific log by name.
If one output is request, it returns the data. If two are requested, it returns the time and the data. Returns empty if the logs don't exist (use the 'contains' function to test for this).
Example:
x = get_log(log, "signal1");
t, x = get_log(log, "signal1");
NOTE that |t| will be ns-by-1 and x will be ns-by-nx, where nx is the number of elements in a single sample.
ODEHybrid.initialize!
— MethodClears out all fields in a Time Series Logger (in-place)
ODEHybrid.interp1
— FunctionSimplified Interp1 (Ported from MATLAB), which returns interpolated values of a 1-D function at query points using linear interpolation.
Inputs: x: Vector containing sample points (in a tuple, pass in (x,) if needs be) v: Vector containing values corresponding to x xq: Vector containing coordinates of the query points options: (Optional) Interpolation method
Outputs: interp_vals: Interpolated values at query points
ODEHybrid.interpd
— MethodAn interpolation function for time series which may have multiple values at a single timestep (e.g., discrete updates). It works like MATLAB's built-in interp1.
xi = interpd(t, x, ti); xi = interpd(t, x, ti, mode);
Inputs:
t Times (n-by-1)
x States (n-by-m)
ti Output times (p-by-1)
mode '-' for left-most values, '+' for right-most values (default)
(NOTE that unlike the MATLAB version, no additional arguments can be passed in yet (e.g., 'nearest'))
Outputs:
xi States corresponding to output times (p-by-m)
Example:
# Create some times and states. Note that there are two states at t=3.
t = [2.76, 2.91, 3, 3, 3.12];
x = [0.2, 0.3, 0.4, 1.1, 1.2];
# Create the desired output times.
ti = (2.8:0.1:3.1);
# Interpolate (t, x) linearly at ti, keeping the right-most values.
xi = interpd(t, x, ti)
# Interpolate (t, x) linearly at ti, keeping the left-most values.
xi = interpd(t, x, ti, '-')
See also: examples_odehybrid.
Online doc: http://www.anuncommonlab.com/doc/odehybrid/interpd.html
Copyright 2014 An Uncommon Lab
ODEHybrid.odehybrid
— FunctionHybrid continuous and discrete propagation.
This function propagates an ordinary differential equation along with a discrete update equation in a manner similar to MATLAB's ode45 function (which only propagates an ordinary differential equation). This is useful for implementing discrete-time controllers or simulating processes that are updated at discrete intervals.
A large number of examples can be found in examples_odehybrid or by entering the following at the command line:
home = fileparts(which('examples_odehybrid'));
web(fullfile(home, 'html', 'examples_odehybrid.html'));
Interfaces:
t, yc, td, yd = odehybrid(solver, ode, de,
dt, ts, yc0, yd0);
t, yc1m, td, yd1n = odehybrid(solver, ode, de,
dt, ts, [yc1, yc2m] [yd1, yd2n]);
t, ..., td, ..., te, yc1m, yd1n, ie = odehybrid(solver, ode, de,
dt, ts, yc0, yd0);
sol = odehybrid(solver, ode, [de1, de2, ...], [dt1, dt2, ...],
ts, yc0, yd0);
sol = odehybrid(..., [options], [log]);
sol = odehybrid(...);
Inputs:
solver Continuous-time propagator to use, e.g. ODEHybrid.rkadapt
ode Ordinary differential equation to use with solver. The interface
should be fun(t, xc1, xc2, ..., xcm, xd1, xd2, ..., xdn) where
xc1, xc2, ..., xcm are the continuous states and xd1, xd2, ...,
xdn are the discrete states. It should return the derivatives of
the continuous states (m outputs).
de Discrete update equation(s) (either a function or cell
array of functions) with the same inputs as ode but
outputing the updated continuous and discrete states (n+m
outputs).
dt Time step(s) of discrete update equation(s). If de is an
array of multiple functions, this should be an array of the same
size.
ts Time span, [t_start, t_end]
yc0 Array of initial continuous states
yd0 Array of initial discrete states
options (Optional) options structure from odeset
log (Optional) TimeSeriesLogger for logging in the ode and de. If a
log is passed in, both ode and de *must* be able to accomodate
having a log input or not as the final argument. E.g., |ode|
will be called as: ode(t, xc1, ..., xd1, ..., xdn) and
ode(t, xc1, ..., xd1, ..., xdn, log).
Outputs:
t Times corresponding to continuous states (nc-by-1)
yc1m Continuous state outputs (m outputs, each nc-by-1)
td Times corresponding to discrete state updates (nd-by-1)
yd1n Discrete state outputs (n outputs, each nd-by-1)
te Times corresponding to events (ne-by-1)
yce1m Continuous states at events (m outputs, each ne-by-1)
yde1n Discrete states at events (n outputs, each ne-by-1)
ie Index of triggering event. See documentation in odeset for more
on events.
sol If a single output is requested, it will be a structure with
fields for each of the individual outputs. The various
continuous states will be grouped in 'yc', the discrete into
'yd', the continuous states at events into 'yce', and the
discrete states at events into 'yde'.
(NOTE that events are not currently supported!)
Example:
# This is a quick example of simulating an unstable continuous system with
# a stabilizing discrete-time controller.
ode = (t, x, u) -> [0 1; 2 0] * x + [0; 1] * u; # Differential equation
de = (t, x, u) -> (x, -[8 4] * x); # Discrete update equation
dt = 0.1; # Discrete eq. time step (float or float array)
ts = [0 5]; # From 0 to 5s
x0 = [1.0; 0.0]; # Initial continuous state (floats, not Ints)
u0 = 0; # Initial discrete state
t, x, tu, u = odehybrid(ODEHybrid.rkadapt, ode, de, dt, ts, x0, u0); # Simulate!
plot(t, x, xlabel = "Time", title = "Example", label = ["x₁" "x₂"])
scatter!(tu, u, label = "u")
See also: examples_odehybrid, ode45, rkadapt, rkfixed.
Online doc: http://www.anuncommonlab.com/doc/odehybrid/odehybrid.html
Copyright 2014 An Uncommon Lab
ODEHybrid.rk4
— MethodRunge-Kutta 4th order integration method.
Implements numerical propagation of an ordinary differential equation from some initial value over the desired range. This function is similar to MATLAB's variable-step ODE propagators (e.g., ode45), but uses a fixed step method. This is useful either when one knows an appropriate step size or when a process is interrupted frequently (ode45 and the similar functions in MATLAB will always make at least a certain number of steps between ts(1) and ts(2), which may be very many more than are necessary).
t, x = rk4(ode, ts, x0, options);
Inputs:
ode Ordinary differential equation function
ts Time span, [t_start, t_end]
x0 Initial state (column vector)
options One can specify an options structure instead of dt
so that this function is compatible with ode45 and its ilk. The
only valid fields are MaxStep (the time step) and OutputFcn
(if dt is passed in rather than options, an empty ODE_SET is initialized
with the provided time step dt)
Outputs:
t Time history
x State history, with each row containing the state corresponding to
the time in the same row of t.
Example:
# Simulate an undamped harmonic oscillator for 10s with a 0.1s time
# step, starting from an initial state of [1; 0]:
ode = (t, x) -> [-x[2]; x[1]]
t, x = rk4(ode, [0 10], [1; 0], 0.1);
plot(t, x);
See also: odehybrid, ode45, odeset, rkfixed, rkadapt.
Online doc: http://www.anuncommonlab.com/doc/odehybrid/rk4.html
Copyright 2014 An Uncommon Lab
ODEHybrid.rkadapt
— MethodRunge-Kutta adaptive-step integration using the specified weights, nodes, and Runge-Kutta matrix (or Dormand-Prince 5(4) by default).
Implements numerical propagation of an ordinary differential equation from some initial value over the desired range. This function is similar to MATLAB's variable-step ODE propagators (e.g., ode45), but is generalized to any Butcher tableau. Unlike ode45, however, there's no inherent minimum number of steps that this function will take; if it can finish integration with a single step, it will, whereas ode45 will always take at least a certain number of steps, no matter how small the time interval and how smooth the dynamics. Therefore, this function is commonly more useful with odehybrid than ode45.
This function is generic for all adaptive-step Runge-Kutta methods. That is, any adaptive-step Runge-Kutta propagator can be created by passing the weightes, nodes, and Runge-Kutta matrix (together, the Butcher tableau) into this function. It will correctly take advantage of the first-same-as-last property, eliminating one function evaluation per step when the table has this property. See the examples below.
t, x = rk4adapt(ode, ts, x0);
t, x = rk4adapt(ode, ts, x0, options);
t, x = rk4adapt(ode, ts, x0, a, b, c p);
t, x = rk4adapt(ode, ts, x0, options, a, b, c, p);
Inputs:
ode Ordinary differential equation function (s.t. ẋ = ode(t, x);)
ts Time span, [t_start, t_end]
x0 Initial state (column vector)
options Options structure from odeset. This function uses the
InitialStep, MaxStep, OutputFcn, RelTol, and AbsTol fields only.
a Runge-Kutta matrix (s-by-s)
b Weights (2-by-s), the top row containing the weights for the
state update
c Nodes (s-by-1)
p Lowest order method represented by b (e.g., for
Runge-Kutta-Fehlberg, this is 4).
(NOTE that if a/b/c/p are not specified, this defaults to Dormand-Prince)
Outputs:
t Time history
x State history, with each row containing the state corresponding to
the time in the same row of t.
Example:
# Simulate an undamped harmonic oscillator for 10s starting from an
# initial state of [1; 0] using the default a, b, c, and p
# (Dormand-Prince 5(4)):
t, x = rkadapt((t,x) -> [-x[2]; x[1]], [0 10], [1; 0]);
plot(t, x);
# Add options.
options = ODE_SET(MaxStep = 0.1, InitialStep = 0.05)
t, x = rkadapt( (t, x) -> [-x[2]; x[1]], [0 10], [1; 0], options)
plot(t, x);
# Now simulate with Bogacki-Shampine 3(2).
a = [ 0 0 0 0;
1/2 0 0 0;
0 3/4 0 0;
2/9 1/3 4/9 0];
b = [ 2/9 1/3 4/9 0; # For the update
7/24 1/4 1/3 1/8]; # For the error prediction
c = [ 0 1/2 3/4 1];
p = 2; # The lower of the two orders
t, x = rkadapt( (t,x) -> [-x[2]; x[1]], [0 10], [1; 0], a, b, c, p);
plot(t, x);
See "Runge-Kutta methods" on Wikipedia for discussion of the Butcher tableau (a, b, and c).
Reference:
Dormand, J. R., and P. J. Prince. "A family of embedded Ruge-Kutta formulae." Journal of Computational and Applied Mathematics 6.1 (1980): 19-26. Web. 22 Apr. 2014.
See also: odehybrid, ode45, odeset, rkfixed.
Online doc: http://www.anuncommonlab.com/doc/odehybrid/rkadapt.html
Copyright 2014 An Uncommon Lab
ODEHybrid.rkfixed
— MethodRunge-Kutta fixed-step integration using the specified weights, nodes, and Runge-Kutta matrix (or the Runge-Kutta 4th order "3/8" method by default).
Implements numerical propagation of an ordinary differential equation from some initial value over the desired range. This function is similar to MATLAB's variable-step ODE propagators (e.g., ode45), but uses a fixed step method. This is useful either when one knows an appropriate step size or when a process is interrupted frequently (ode45 and the similar functions in MATLAB will always make at least a certain number of steps between ts(1) and ts(2), which may be very many more than are necessary).
This function is generic for all fixed-step Runge-Kutta methods. That is, any fixed-step Runge-Kutta propagator can be created by passing the weightes, nodes, and Runge-Kutta matrix (together, the Butcher tableau) into this function. See the example below.
t, x = rkfixed(ode, ts, x0, dt);
t, x = rkfixed(ode, ts, x0, dt, a, b, c);
t, x = rkfixed(ode, ts, x0, options);
t, x = rkfixed(ode, ts, x0, options, a, b, c);
Inputs:
ode Ordinary differential equation function
ts Time span, [t_start, t_end]
x0 Initial state (column vector)
options One can specify an options structure instead of dt
so that this function is compatible with ode45 and its ilk. The
only valid fields are MaxStep (the time step) and OutputFcn
a Runge-Kutta matrix
b Weights
c Nodes
Outputs:
t Time history
x State history, with each row containing the state corresponding to
the time in the same row of t.
Example:
# Simulate an undamped harmonic oscillator for 10s with a 0.1s time
# step, starting from an initial state of [1; 0] using RK 4th order
# integration (via the Butcher tableau specified by a, b, and c). This
# is exactly the same as the rk4 function
a = [ 0 0 0 0;
0.5 0 0 0;
0 0.5 0 0;
0 0 1 0];
b = [ 1 2 2 1]/6;
c = [ 0 0.5 0.5 1];
t, x = rkfixed( (t,x) -> [-x[2]; x[1]], [0 10], [1; 0], 0.1, a, b, c);
plot(t, x);
See "Runge-Kutta methods" on Wikipedia for discussion of the Butcher tableau (a, b, and c).
See also: odehybrid, ode45, odeset, rk4.
Online doc: http://www.anuncommonlab.com/doc/odehybrid/rkfixed.html
Copyright 2014 An Uncommon Lab
ODEHybrid.state_to_vector
— MethodBuild a state vector by recursively moving through the elements of a more complex state (matrix, cell array, or struct), saving the numerical values along the way in a vector. This is used in odehybrid to convert a complex series of states into a vector to be used for continuous state updates. The reverse function is vectortostate.
vector = statetovector(state);
Inputs:
state: A numeric, array, or struct array type, the contents of
which consist of numeric, array, or struct array types, etc.
(Structs MUST be mutable)
Outputs:
vector: A vector containing the numerical values from the state (doubles)
Example:
mutable struct TestStruc
a
bcd
end
Base.:(==)(a::TestStruc, b::TestStruc) = (Base.:(==)(a.a, b.a) && Base.:(==)(a.bcd, b.bcd)) # Custom
x = [ [1 3; 4 2], TestStruc([5; 6], [7:9; 10]), 3.14 ]
v = state_to_vector(x)
x2 = vector_to_state(v, x)
x .== x2[1]
See also: vectortostate.
Online doc: http://www.anuncommonlab.com/doc/odehybrid/statetovector.html
Copyright 2014 An Uncommon Lab
ODEHybrid.tsl_plot
— MethodPlot all of the signals, grouped appropriately into figures. A custom x label can be added to figures as well
ODEHybrid.vec_to_mat
— MethodHelper function that converts a N-element vector with each element of length M into an N x M matrix
ODEHybrid.vector_to_state
— MethodBuild a state from a state vector by recursively moving through the elements of a more complex example state (matrix, cell array, or struct), placing the relevant numerical values from the vector into the appropriate places. This is used in odehybrid to convert a state vector into the original complex setup of states. The reverse function is statetovector.
state = vectortostate(vector, state);
Inputs:
vector A vector containing the numerical values from the object
state A numeric, cell array, or struct array type representing the
sizes and types for the values in vector
Outputs:
state A numeric, cell array, or struct array type
count (Internal use only)
Example:
mutable struct TestStruc
a
bcd
end
Base.:(==)(a::TestStruc, b::TestStruc) = (Base.:(==)(a.a, b.a) && Base.:(==)(a.bcd, b.bcd)) # Custom
x = [ [1 3; 4 2], TestStruc([5; 6], [7:9; 10]), 3.14 ]
v = state_to_vector(x)
x2 = vector_to_state(v, x)
x .== x2[1]
See also: statetovector.
Online doc: http://www.anuncommonlab.com/doc/odehybrid/vectortostate.html
Copyright 2014 An Uncommon Lab