dysts API Reference¶
The API reference for the dysts .. _repository: https://github.com/williamgilpin/dysts
Base Classes¶
Dynamical systems in Python
(M, T, D) or (T, D) convention for outputs
Requirements: + numpy + scipy + sdeint (for integration with noise) + numba (optional, for faster integration)
- class dysts.base.BaseDyn(**entries)¶
A base class for dynamical systems
- name¶
The name of the system
- Type:
str
- params¶
The parameters of the system.
- Type:
dict
- random_state¶
The seed for the random number generator. Defaults to None
- Type:
int
- Development:
Add a function to look up additional metadata, if requested
- static bound_trajectory(traj)¶
Bound a trajectory within a periodic domain
- load_trajectory(subsets='train', granularity='fine', return_times=False, standardize=False, noise=False)¶
Load a precomputed trajectory for the dynamical system
- Parameters:
subsets ("train" | "test") – Which dataset (initial conditions) to load
granularity ("course" | "fine") – Whether to load fine or coarsely-spaced samples
noise (bool) – Whether to include stochastic forcing
standardize (bool) – Standardize the output time series.
return_times (bool) – Whether to return the timepoints at which the solution was computed
- Returns:
A T x D trajectory tpts, sol (ndarray): T x 1 timepoint array, and T x D trajectory
- Return type:
sol (ndarray)
- make_trajectory(*args, **kwargs)¶
Make a trajectory for the dynamical system
- sample(*args, **kwargs)¶
Sample a trajectory for the dynamical system via numerical integration
- update_params()¶
Update all instance attributes to match the values stored in the params field
- class dysts.base.DynMap(**kwargs)¶
A dynamical system base class, which loads and assigns parameter values from a file
- Parameters:
params (list) – parameter values for the differential equations
kwargs (dict) – A dictionary of keyword arguments passed to the base dynamical model class
- make_trajectory(n, inverse=False, return_times=False, standardize=False, **kwargs)¶
Generate a fixed-length trajectory with default timestep, parameters, and initial condition(s)
- Parameters:
n (int) – the length of each trajectory
inverse (bool) – whether to reverse a trajectory
standardize (bool) – Standardize the output time series.
return_times (bool) – Whether to return the timepoints at which the solution was computed
- rhs(X)¶
The right hand side of a dynamical map
- rhs_inv(Xp)¶
The inverse of the right hand side of a dynamical map
- class dysts.base.DynSys(**kwargs)¶
A continuous dynamical system base class, which loads and assigns parameter values from a file
- kwargs¶
A dictionary of keyword arguments passed to the base dynamical model class
- Type:
dict
- make_trajectory(n, method='Radau', resample=True, pts_per_period=100, return_times=False, standardize=False, postprocess=True, noise=0.0)¶
Generate a fixed-length trajectory with default timestep, parameters, and initial conditions
- Parameters:
n (int) – the total number of trajectory points
method (str) – the integration method
resample (bool) – whether to resample trajectories to have matching dominant Fourier components
pts_per_period (int) – if resampling, the number of points per period
standardize (bool) – Standardize the output time series.
return_times (bool) – Whether to return the timepoints at which the solution was computed
postprocess (bool) – Whether to apply coordinate conversions and other domain-specific rescalings to the integration coordinates
noise (float) – The amount of stochasticity in the integrated dynamics. This would correspond to Brownian motion in the absence of any forcing.
- Returns:
A T x D trajectory tpts, sol (ndarray): T x 1 timepoint array, and T x D trajectory
- Return type:
sol (ndarray)
- rhs(X, t)¶
The right hand side of a dynamical equation
- class dysts.base.DynSysDelay(**kwargs)¶
A delayed differential equation object. Defaults to using Euler integration scheme The delay timescale is assumed to be the “tau” field. The embedding dimension is set by default to ten, but delay equations are infinite dimensional. Uses a double-ended queue for memory efficiency
- kwargs¶
A dictionary of keyword arguments passed to the dynamical system parent class
- Type:
dict
- make_trajectory(n, d=10, method='Euler', noise=0.0, resample=False, pts_per_period=100, standardize=False, return_times=False, postprocess=True)¶
Generate a fixed-length trajectory with default timestep, parameters, and initial conditions.
- Parameters:
n (int) – the total number of trajectory points
d (int) – the number of embedding dimensions to return
method (str) – Not used. Currently Euler is the only option here
noise (float) – The amplitude of brownian forcing
resample (bool) – whether to resample trajectories to have matching dominant Fourier components
pts_per_period (int) – if resampling, the number of points per period
standardize (bool) – Standardize the output time series.
return_times (bool) – Whether to return the timepoints at which the solution was computed
- rhs(X, t)¶
The right hand side of a dynamical equation
- dysts.base.get_attractor_list(model_type='continuous')¶
Returns the names of all models in the package
- Parameters:
model_type (str) – “continuous” (default) or “discrete”
- Returns:
The names of all attractors in database
- Return type:
attractor_list (list of str)
- dysts.base.make_trajectory_ensemble(n, subset=None, use_multiprocessing=False, random_state=None, **kwargs)¶
Integrate multiple dynamical systems with identical settings
- Parameters:
n (int) – The number of timepoints to integrate
subset (list) – A list of system names. Defaults to all systems
use_multiprocessing (bool) – Not yet implemented.
random_state (int) – The random seed to use for the ensemble
kwargs (dict) – Integration options passed to each system’s make_trajectory() method
- Returns:
A dictionary containing trajectories for each system
- Return type:
all_sols (dict)
Datasets¶
- class dysts.datasets.TimeSeriesDataset(datapath=None, data=None)¶
A data structure for handling time series. Reads and writes from JSON
Assumes that all time series have the same number of timepoints. In the case of multivariate series with different numbers of features / dimensions, the missing dimensionality is padded
- Parameters:
data (dict, optional) – A time series dataset, consisting of a dictionary with index given by top-level names of different time series.
datapath (string, optional) – The path to a JSON dataset. Defaults to creating an empty dataset object
- rank¶
The original length of each time series
- Type:
list
- is_multivariate¶
Whether the time series dataset is multivariate
- Type:
bool
- names¶
The names of the systems in the dataset
- Type:
list
- get_rowvalues(key)¶
Retrieve the value of a field from each row
- to_array(standardize=False)¶
Return the values of a time series as a matrix of shape (B, T, D) or (B, T)
- to_pandas(standardize=True, filling='constant')¶
Convert a data dictionary to a DataFrame Draws values from self.data (dict), a dictionary with index given by top-level keys
- Parameters:
standardize (bool) – whether to standardize each dataset before flattening
filling (str) – For multivariate time series of varying demensionality, values to pad in empty columns.
- Returns:
a 2D dataset consisting of indexed time series
- Return type:
data_pd (pd.DataFrame)
- trim_series(start_val, stop_val)¶
Trims a time series in-place. This is not reversible without re-initializing
- dysts.datasets.convert_json_to_gzip(fpath, encoding='utf-8')¶
Convert a json file to a gzip file in a format that can be easily read by the dysts package. By default, the gzip file will be saved with the same name and in the same directory as the json file, but with a “.gz” extension.
- Parameters:
fpath (str) – Path to the json file to be converted
encoding (str) – Encoding to use when writing the gzip file
- Returns:
None
- dysts.datasets.load_dataset(subsets='train', univariate=True, granularity='fine', data_format='object', noise=False, split_fraction=0.8333333333333334, **kwargs)¶
Load dynamics from continuous dynamical systems.
- Parameters:
subsets ("train" | "train_val" | "test" | "test_val" | "train_all" | "test_all") – Which dataset to draw. Train and train val correspond to the same time series split 5/6 of the way through, while “test” and “test_val” both represent a trajectory emanating from a different initial condition, split 5/6 of the way though. “train_all” and “test_all” return full trajectories without splits
univariate (bool) – Whether to use one coordinate, or all for each system.
granularity ("course" | "fine") – Whether to use fine or coarsely-spaced samples
data_format ("object" | "numpy" | "pandas") – The format to return
noise (bool) – Whether to include stochastic forcing
split_fraction (float) – The fraction of the time series to hold out as a test/val partition
kwargs (dict) – keyword arguments passed to the data formatter
- Returns:
A collection of time series datasets
- Return type:
dataset (TimeSeriesDataset)
Utilities¶
Helper utilities for working with time series arrays. This module is intended to have no dependencies on the rest of the package.
- Author: Jacob Reinhold (
)
Created on: June 10, 2020
- Copyright: 2020 Jacob Reinhold (
)
- class dysts.utils.ComputationHolder(func=None, *args, timeout=10, **kwargs)¶
A wrapper class to force a computation to stop after a timeout.
- Parameters
func (callable): the function to run args (tuple): the arguments to pass to the function kwargs (dict): the keyword arguments to pass to the function timeout (int): the timeout in seconds. If None is passed, the computation
will run indefinitely until it finishes.
- Example
>>> def my_func(): ... while True: ... print("hello") ... time.sleep(8) >>> ch = ComputationHolder(my_func, timeout=3) >>> ch.run() hello hello hello None
- dysts.utils.cartesian_to_polar(x, y)¶
Convert 2D cartesian coordinates to polar coordinates
- dysts.utils.dict_loudassign(d, key, val)¶
Assign a key val pair in a dict, and print the result
- dysts.utils.find_characteristic_timescale(y, k=1, window=True)¶
Find the k leading characteristic timescales in a time series using the power spectrum..
- dysts.utils.find_psd(y, window=True)¶
Find the power spectrum of a signal
- dysts.utils.find_significant_frequencies(sig, window=True, fs=1, n_samples=100, significance_threshold=0.95, surrogate_method='rp', show=False, return_amplitudes=False)¶
Find power spectral frequencies that are significant in a signal, by comparing the appearance of a peak with its appearance in randomly-shuffled surrogates.
If no significant freqencies are detected, the significance floor is lowered
- Parameters:
window (bool) – Whether to window the signal before taking the FFT
thresh (float) – The number of standard deviations above mean to be significant
fs (int) – the sampling frequency
n_samples (int) – the number of surrogates to create
show (bool) – whether to show the psd of the signal and the surrogate
- Returns:
The frequencies overrated in the dataset amps (ndarray): the amplitudes of the PSD at the identified frequencies
- Return type:
freqs (ndarray)
- dysts.utils.find_slope(x, y)¶
Given two vectors or arrays, compute the best fit slope using an analytic formula. For arrays is computed along the last axis.
- Parameters:
x (ndarray) – (N,) or (M, N)
y (ndarray) – (N,) or (M, N)
- Returns:
the values of the slope for each of the last dimensions
- Return type:
b (ndarray)
- dysts.utils.freq_from_autocorr(sig, fs=1)¶
Estimate frequency using autocorrelation
- Parameters:
sig (ndarray) – A univariate signal
fs (int) – The sampling frequency
- Returns:
The dominant frequency
- Return type:
out (float)
References
Modified from the following https://gist.github.com/endolith/255291
- dysts.utils.freq_from_fft(sig, fs=1)¶
Estimate frequency of a signal from the peak of the power spectrum
- Parameters:
sig (ndarray) – A univariate signal
fs (int) – The sampling frequency
- Returns:
The dominant frequency
- Return type:
out (float)
References
Modified from the following https://gist.github.com/endolith/255291
- dysts.utils.generate_ic_ensemble(model, tpts0, n_samples, frac_perturb_param=0.1, frac_transient=0.1, ic_range=None, random_state=0)¶
Generate an ensemble of trajectories with random initial conditions, labelled by different initial conditions
- Parameters:
model (callable_) – function defining the numerical derivative
tpts (ndarray) – the timesteps over which to run the simulation
n_samples (int) – the number of different initial conditons
frac_perturb_param (float) – the amount to perturb the ic by
frac_transient (float) – the fraction of time for the time series to settle onto the attractor
ic_range (list) – a starting value for the initial conditions
random_state (int) – the seed for the random number generator
- Returns:
A set of initial conditions
- Return type:
all_samples (array)
- dysts.utils.integrate_dyn(f, ic, tvals, noise=0, dtval=None, **kwargs)¶
Given the RHS of a dynamical system, integrate the system noise > 0 requires the Python library sdeint (assumes Brownian noise)
- Parameters:
f (callable) – The right hand side of a system of ODEs.
ic (ndarray) – the initial conditions
noise_amp (float or iterable) – The amplitude of the Langevin forcing term. If a vector is passed, this will be different for each dynamical variable
dtval (float) – The starting integration timestep. This will be the exact timestep for fixed-step integrators, or stochastic integration.
kwargs (dict) – Arguments passed to scipy.integrate.solve_ivp.
- Returns:
The integrated trajectory
- Return type:
sol (ndarray)
- dysts.utils.jac_fd(func0, y0, eps=0.001, m=1, method='central', verbose=False)¶
Calculate numerical jacobian of a function with respect to a reference value
- Parameters:
func (callable) – a vector-valued function
y0 (ndarray) – a point around which to take the gradient
eps (float) – the step size for the finite difference calculation
- Returns:
a numerical estimate of the Jacobian about that point
- Return type:
jac (ndarray)
- dysts.utils.make_epsilon_ball(pt, n, eps=1e-05, random_state=None)¶
Uniformly sample a fixed-radius ball of points around a given point via using Muller’s method
- Parameters:
pt (ndarray) – The center of the sampling
n (int) – The number of points to sample
eps (float) – The radius of the ball
random_state (int) – Initialize the random number generator
- Returns:
The set of randomly-sampled points
- Return type:
out (ndarray)
- dysts.utils.make_surrogate(data, method='rp')¶
- Parameters:
data (ndarray) – A one-dimensional time series
method (str) – “rs” or rp”
- Returns:
A single random surrogate time series
- Return type:
surr_data (ndarray)
- dysts.utils.pad_axis(arr, d, axis=-1, padding=0)¶
Pad axis of arr with a constant padding to a desired shape
- dysts.utils.pad_to_shape(arr, target_shape)¶
Given an array, and a target shape, pad the dimensions in order to reach the desired shape Currently, if the rank of the array is lower than the target shape, singleton dimensions are appended to the rank
- Parameters:
arr (ndarray) – The array to pad.
target_shape (iterable) – The desired shape.
- Returns:
The padded array,
- Return type:
arr (ndarray)
- dysts.utils.polar_to_cartesian(r, th)¶
Convert polar coordinates to 2D Cartesian coordinates
- dysts.utils.standardize_ts(a, scale=1.0)¶
Standardize an array along dimension -2 For dimensions with zero variance, divide by one instead of zero
- Parameters:
a (ndarray) – a matrix containing a time series or batch of time series with shape (T, D) or (B, T, D)
scale (float) – the number of standard deviations by which to scale
- Returns:
- A standardized time series with the same shape as
the input
- Return type:
ts_scaled (ndarray)
Analysis¶
Functions that act on DynSys or DynMap objects
- dysts.analysis.calculate_lyapunov_exponent(traj1, traj2, dt=1.0)¶
Calculate the lyapunov exponent of two multidimensional trajectories using simple linear regression based on the log-transformed separation of the trajectories.
- Parameters:
traj1 (np.ndarray) – trajectory 1 with shape (n_timesteps, n_dimensions)
traj2 (np.ndarray) – trajectory 2 with shape (n_timesteps, n_dimensions)
dt (float) – time step between timesteps
- Returns:
lyapunov exponent
- Return type:
float
- dysts.analysis.compute_timestep(model, total_length=40000, transient_fraction=0.2, num_iters=20, pts_per_period=1000, visualize=False, return_period=True)¶
Given a dynamical system object, find the integration timestep based on the largest signficant frequency
- Parameters:
model (DynSys) – A dynamical systems object.
total_length (int) – The total trajectory length to use to determine timescales.
transient_fraction (float) – The fraction of a trajectory to discard as transient
num_iters (int) – The number of refinements to the timestep
pts_per_period (int) – The target integration timestep relative to the signal.
visualize (bool) – Whether to plot timestep versus time, in order to identify problems with the procedure
return_period (bool) – Whether to calculate and retunr the dominant timescale in the signal
- Returns
dt (float): The best integration timestep period (float, optional): The dominant timescale in the signal
- dysts.analysis.find_lyapunov_exponents(model, traj_length, pts_per_period=500, tol=1e-08, min_tpts=10, **kwargs)¶
Given a dynamical system, compute its spectrum of Lyapunov exponents. :param model: the right hand side of a differential equation, in format
func(X, t)
- Parameters:
traj_length (int) – the length of each trajectory used to calulate Lyapunov exponents
pts_per_period (int) – the sampling density of the trajectory
kwargs – additional keyword arguments to pass to the model’s make_trajectory method
- Returns:
A list of computed Lyapunov exponents
- Return type:
final_lyap (ndarray)
References
- Christiansen & Rugh (1997). Computing Lyapunov spectra with continuous
Gram-Schmidt orthonormalization
Example
>>> import dysts >>> model = dysts.Lorenz() >>> lyap = dysts.find_lyapunov_exponents(model, 1000, pts_per_period=1000) >>> print(lyap)
- dysts.analysis.get_train_test(eq, n_train=1000, n_test=200, standardize=True, **kwargs)¶
Generate train and test trajectories for a given dynamical system
- Parameters:
eq (dysts.DynSys) – a dynamical system object
n_train (int) – number of points in the training trajectory
n_test (int) – number of points in the test trajectory
standardize (bool) – whether to standardize the trajectories
**kwargs – additional keyword arguments to pass to make_trajectory
- Returns:
- a tuple containing:
- (tuple): a tuple containing:
(ndarray): the timepoints of the training trajectory (ndarray): the training trajectory
- (tuple): a tuple containing:
(ndarray): the timepoints of the test trajectory (ndarray): the test trajectory
- Return type:
(tuple)
- dysts.analysis.kaplan_yorke_dimension(spectrum0)¶
Calculate the Kaplan-Yorke dimension, given a list of Lyapunov exponents
- dysts.analysis.lyapunov_exponent_naive(eq, rtol=0.001, atol=1e-10, n_samples=1000, traj_length=5000, max_walltime=None, **kwargs)¶
Calculate the lyapunov spectrum of the system using a naive method based on the log-transformed separation of the trajectories over time.
- Parameters:
eq (dysts.DynSys) – equation to calculate the lyapunov spectrum of
rtol (float) – relative tolerance for the separation of the trajectories
atol (float) – absolute tolerance for the separation of the trajectories
n_samples (int) – number of initial conditions to sample
traj_length (int) – length of the trajectories to sample. This should be long enough to ensure that most trajectories explore the attractor.
max_walltime (float) – maximum walltime in seconds to spend on the calculation of a given trajectory. If the calculation takes longer than this, the trajectory is discarded and a new one is sampled.
**kwargs – keyword arguments to pass to the sample / make_trajectory method of the dynamical equation
- Returns:
largest lyapunov exponent
- Return type:
float
Example
>>> import dysts >>> eq = dysts.Lorenz() >>> max_lyap = dysts.lyapunov_exponent_naive(eq)
- dysts.analysis.mse_mv(traj)¶
Generate an estimate of the multivariate multiscale entropy. The current version computes the entropy separately for each channel and then averages. It therefore represents an upper-bound on the true multivariate multiscale entropy
- Parameters:
traj (ndarray) – a trajectory of shape (n_timesteps, n_channels)
- Returns:
the multivariate multiscale entropy
- Return type:
mmse (float)
- dysts.analysis.sample_initial_conditions(model, points_to_sample, traj_length=1000, pts_per_period=30)¶
Generate a random sample of initial conditions from a dynamical system
- Parameters:
model (callable) – the right hand side of a differential equation, in format func(X, t)
points_to_sample (int) – the number of random initial conditions to sample
traj_length (int) – the total length of the reference trajectory from which points are drawn
pts_per_period (int) – the sampling density of the trajectory
- Returns:
The points with shape (points_to_sample, d)
- Return type:
sample_points (ndarray)