pykoopman.common package¶
Submodules¶
pykoopman.common.cqgle module¶
Module for cubic-quintic Ginzburg-Landau equation.
- class pykoopman.common.cqgle.cqgle(n, x, dt, tau=0.08, kappa=0, beta=0.66, nu=-0.1, sigma=-0.1, gamma=-0.1, L=6.283185307179586)[source]¶
Bases:
object
Cubic-quintic Ginzburg-Landau equation solver.
Solves the equation: i*u_t + (0.5 - i * tau) u_{xx} - i * kappa u_{xxxx} + (1-i * beta)|u|^2 u + (nu - i * sigma)|u|^4 u - i * gamma u = 0
Solves the periodic boundary conditions PDE using spectral methods.
- Attributes:
n_states (int) – Number of states.
x (numpy.ndarray) – x-coordinates.
dt (float) – Time step.
tau (float) – Parameter tau.
kappa (float) – Parameter kappa.
beta (float) – Parameter beta.
nu (float) – Parameter nu.
sigma (float) – Parameter sigma.
gamma (float) – Parameter gamma.
k (numpy.ndarray) – Wave numbers.
dk (float) – Wave number spacing.
- collect_data_continuous(x0)[source]¶
collect training data pairs - continuous sense.
given x0, with shape (n_dim, n_traj), the function returns dx/dt with shape (n_dim, n_traj)
pykoopman.common.examples module¶
module for example dynamics data
- pykoopman.common.examples.drss(n=2, p=2, m=2, p_int_first=0.1, p_int_others=0.01, p_repeat=0.05, p_complex=0.5)[source]¶
Create a discrete-time, random, stable, linear state space model.
- Parameters:
n (int, optional) – Number of states. Default is 2.
p (int, optional) – Number of control inputs. Default is 2.
m (int, optional) – Number of output measurements. If m=0, C becomes the identity matrix, so that y=x. Default is 2.
p_int_first (float, optional) – Probability of an integrator as the first pole. Default is 0.1.
p_int_others (float, optional) – Probability of other integrators beyond the first. Default is 0.01.
p_repeat (float, optional) – Probability of repeated roots. Default is 0.05.
p_complex (float, optional) – Probability of complex roots. Default is 0.5.
- Returns:
A tuple containing the state transition matrix (A), control matrix (B), and measurement matrix (C).
A (numpy.ndarray): State transition matrix of shape (n, n). B (numpy.ndarray): Control matrix of shape (n, p). C (numpy.ndarray): Measurement matrix of shape (m, n). If m = 0, C is the
identity matrix.
- Return type:
Tuple[numpy.ndarray, numpy.ndarray, numpy.ndarray]
- pykoopman.common.examples.advance_linear_system(x0, u, n, A=None, B=None, C=None)[source]¶
Simulate the linear system dynamics for a given number of steps.
- Parameters:
x0 (numpy.ndarray) – Initial state vector of shape (n,).
u (numpy.ndarray) – Control input array of shape (p,) or (p, n-1). If 1-dimensional, it will be converted to a row vector.
n (int) – Number of steps to simulate.
A (numpy.ndarray, optional) – State transition matrix of shape (n, n). If not provided, it defaults to None.
B (numpy.ndarray, optional) – Control matrix of shape (n, p). If not provided, it defaults to None.
C (numpy.ndarray, optional) – Measurement matrix of shape (m, n). If not provided, it defaults to None.
- Returns:
A tuple containing the state trajectory (x) and the output trajectory (y).
x (numpy.ndarray): State trajectory of shape (n, len(x0)). y (numpy.ndarray): Output trajectory of shape (n, C.shape[0]).
- Return type:
Tuple[numpy.ndarray, numpy.ndarray]
- pykoopman.common.examples.vdp_osc(t, x, u)[source]¶
Compute the dynamics of the Van der Pol oscillator.
- Parameters:
t (float) – Time.
x (numpy.ndarray) – State vector of shape (2,).
u (float) – Control input.
- Returns:
Updated state vector of shape (2,).
- Return type:
numpy.ndarray
- pykoopman.common.examples.rk4(t, x, u, _dt=0.01, func=<function vdp_osc>)[source]¶
Perform a 4th order Runge-Kutta integration.
- Parameters:
t (float) – Time.
x (numpy.ndarray) – State vector of shape (2,).
u (float) – Control input.
_dt (float, optional) – Time step. Defaults to 0.01.
func (function, optional) – Function defining the dynamics. Defaults to vdp_osc.
- Returns:
Updated state vector of shape (2,).
- Return type:
numpy.ndarray
- pykoopman.common.examples.square_wave(step)[source]¶
Generate a square wave with a period of 60 time steps.
- Parameters:
step (int) – Current time step.
- Returns:
Square wave value at the given time step.
- Return type:
float
- pykoopman.common.examples.sine_wave(step)[source]¶
Generate a sine wave with a period of 60 time steps.
- Parameters:
step (int) – Current time step.
- Returns:
Sine wave value at the given time step.
- Return type:
float
- pykoopman.common.examples.lorenz(x, t, sigma=10, beta=2.6666666666666665, rho=28)[source]¶
Compute the derivative of the Lorenz system at a given state.
- Parameters:
x (list) – Current state of the Lorenz system [x, y, z].
t (float) – Current time.
sigma (float, optional) – Parameter sigma. Default is 10.
beta (float, optional) – Parameter beta. Default is 8/3.
rho (float, optional) – Parameter rho. Default is 28.
- Returns:
Derivative of the Lorenz system [dx/dt, dy/dt, dz/dt].
- Return type:
list
- pykoopman.common.examples.rev_dvdp(t, x, u=0, dt=0.1)[source]¶
Reverse dynamics of the Van der Pol oscillator.
- Parameters:
t (float) – Time.
x (numpy.ndarray) – Current state of the system [x1, x2].
u (float, optional) – Input. Default is 0.
dt (float, optional) – Time step. Default is 0.1.
- Returns:
Updated state of the system [x1’, x2’].
- Return type:
numpy.ndarray
- class pykoopman.common.examples.Linear2Ddynamics[source]¶
Bases:
object
Initializes a Linear2Ddynamics object.
- linear_map(x)[source]¶
Applies the linear mapping to the input state.
- Parameters:
x (numpy.ndarray) – Input state.
- Returns:
Resulting mapped state.
- Return type:
numpy.ndarray
- class pykoopman.common.examples.torus_dynamics(n_states=128, sparsity=5, freq_max=15, noisemag=0.0)[source]¶
Bases:
object
Initializes a torus_dynamics object.
- Parameters:
n_states (int, optional) – Number of states. Default is 128.
sparsity (int, optional) – Degree of sparsity. Default is 5.
freq_max (int, optional) – Maximum frequency. Default is 15.
noisemag (float, optional) – Magnitude of noise. Default is 0.0.
- advance(n_samples, dt=1)[source]¶
Advances the continuous-time dynamics without control.
- Parameters:
n_samples (int) – Number of samples to advance.
dt (float, optional) – Time step. Default is 1.
- advance_discrete_time(n_samples, dt, u=None)[source]¶
Advances the discrete-time dynamics with or without control.
- Parameters:
n_samples (int) – Number of samples to advance.
dt (float) – Time step.
u (array-like, optional) – Control input. Default is None.
- set_control_matrix_physical(B)[source]¶
Sets the control matrix in physical space.
- Parameters:
B (array-like) – Control matrix in physical space.
- set_control_matrix_fourier(Bhat)[source]¶
Sets the control matrix in Fourier space.
- Parameters:
Bhat (array-like) – Control matrix in Fourier space.
- set_point_actuator(position=None)[source]¶
Sets a single point actuator.
- Parameters:
position (array-like, optional) – Position of the actuator. Default is None.
- viz_torus(ax, x)[source]¶
Visualizes the torus dynamics.
- Parameters:
ax – Axes object for plotting.
x (array-like) – Dynamics to be visualized.
- Returns:
Surface plot of the torus dynamics.
- Return type:
surface
- viz_all_modes(modes=None)[source]¶
Visualizes all modes.
- Parameters:
modes (array-like, optional) – Modes to be visualized. Default is None.
- Returns:
Figure object containing the visualizations.
- Return type:
fig
- property modes¶
Returns the modes of the dynamics.
- Returns:
Modes of the dynamics.
- Return type:
modes (array-like)
- property B_effective¶
Returns the effective control matrix.
- Returns:
Effective control matrix.
- Return type:
B_effective (array-like)
- class pykoopman.common.examples.slow_manifold(mu=-0.05, la=-1.0, dt=0.01)[source]¶
Bases:
object
Represents the slow manifold class.
- Parameters:
mu (float, optional) – Parameter mu. Default is -0.05.
la (float, optional) – Parameter la. Default is -1.0.
dt (float, optional) – Time step size. Default is 0.01.
- Attributes:
mu (float) – Parameter mu.
la (float) – Parameter la.
b (float) – Value computed from mu and la.
dt (float) – Time step size.
n_states (int) – Number of states.
- sys(t, x, u)[source]¶
Computes the system dynamics.
- Parameters:
t (float) – Time.
x (array-like) – State.
u (array-like) – Control input.
- Returns:
Computed system dynamics.
- Return type:
array-like
- output(x)[source]¶
Computes the output based on the state.
- Parameters:
x (array-like) – State.
- Returns:
Computed output.
- Return type:
array-like
- simulate(x0, n_int)[source]¶
Simulates the system dynamics.
- Parameters:
x0 (array-like) – Initial state.
n_int (int) – Number of integration steps.
- Returns:
Simulated trajectory.
- Return type:
array-like
- collect_data_continuous(x0)[source]¶
Collects data from continuous-time dynamics.
- Parameters:
x0 (array-like) – Initial state.
- Returns:
Collected data (X, Y).
- Return type:
tuple
- collect_data_discrete(x0, n_int)[source]¶
Collects data from discrete-time dynamics.
- Parameters:
x0 (array-like) – Initial state.
n_int (int) – Number of integration steps.
- Returns:
Collected data (X, Y).
- Return type:
tuple
- class pykoopman.common.examples.forced_duffing(dt, d, alpha, beta)[source]¶
Bases:
object
Initializes the Forced Duffing Oscillator.
- Parameters:
dt (float) – Time step.
d (float) – Damping coefficient.
alpha (float) – Coefficient of x1.
beta (float) – Coefficient of x1^3.
- sys(t, x, u)[source]¶
Defines the system dynamics of the Forced Duffing Oscillator.
- Parameters:
t (float) – Time.
x (array-like) – State vector.
u (array-like) – Control input.
- Returns:
Rate of change of the state vector.
- Return type:
array-like
- simulate(x0, n_int, u)[source]¶
Simulates the Forced Duffing Oscillator.
- Parameters:
x0 (array-like) – Initial state vector.
n_int (int) – Number of time steps.
u (array-like) – Control inputs.
- Returns:
State trajectories.
- Return type:
array-like
- collect_data_continuous(x0, u)[source]¶
Collects continuous data for the Forced Duffing Oscillator.
- Parameters:
x0 (array-like) – Initial state vector.
u (array-like) – Control inputs.
- Returns:
State and output trajectories.
- Return type:
tuple
- collect_data_discrete(x0, n_int, u)[source]¶
Collects discrete-time data for the Forced Duffing Oscillator.
- Parameters:
x0 (array-like) – Initial state vector.
n_int (int) – Number of time steps.
u (array-like) – Control inputs.
- Returns:
State and output trajectories.
- Return type:
tuple
pykoopman.common.ks module¶
module for 1D KS equation
- class pykoopman.common.ks.ks(n, x, nu, dt, M=16)[source]¶
Bases:
object
Solving 1D KS equation
u_t = -u*u_x + u_{xx} + nu*u_{xxxx}
Periodic B.C. between 0 and 2*pi. This PDE is solved using spectral methods.
pykoopman.common.nlse module¶
module for nonlinear schrodinger equation
- class pykoopman.common.nlse.nlse(n, dt, L=6.283185307179586)[source]¶
Bases:
object
nonlinear schrodinger equation
iu_t + 0.5u_xx + u*|u|^2 = 0
periodic B.C. PDE is solved with Spectral methods using FFT
- collect_data_continuous(x0)[source]¶
collect training data pairs - continuous sense.
given x0, with shape (n_dim, n_traj), the function returns dx/dt with shape (n_dim, n_traj)
pykoopman.common.validation module¶
- pykoopman.common.validation.drop_nan_rows(arr, *args)[source]¶
Remove rows in all inputs for which
arr
has_np.nan
entries.- Parameters:
arr (numpy.ndarray) – Array whose rows are checked for nan entries. Any rows containing nans are removed from
arr
and all arguments passed viaargs
.*args (variable length argument list of numpy.ndarray) – Additional arrays from which to remove rows. Each argument should have the same number of rows as
arr
.
- Returns:
arrays – Arrays with nan rows dropped. The first entry corresponds to
arr
and all following entries to*args
.- Return type:
tuple of numpy.ndarray
pykoopman.common.vbe module¶
module for 1D viscous burgers
- class pykoopman.common.vbe.vbe(n, x, dt, nu=0.1, L=6.283185307179586)[source]¶
Bases:
object
1D viscous Burgers equation
u_t = -u*u_x +
u u_{xx}
periodic B.C. PDE is solved using spectral methods
- collect_data_continuous(x0)[source]¶
collect training data pairs - continuous sense.
given x0, with shape (n_dim, n_traj), the function returns dx/dt with shape (n_dim, n_traj)
Module contents¶
- pykoopman.common.drop_nan_rows(arr, *args)[source]¶
Remove rows in all inputs for which
arr
has_np.nan
entries.- Parameters:
arr (numpy.ndarray) – Array whose rows are checked for nan entries. Any rows containing nans are removed from
arr
and all arguments passed viaargs
.*args (variable length argument list of numpy.ndarray) – Additional arrays from which to remove rows. Each argument should have the same number of rows as
arr
.
- Returns:
arrays – Arrays with nan rows dropped. The first entry corresponds to
arr
and all following entries to*args
.- Return type:
tuple of numpy.ndarray
- pykoopman.common.drss(n=2, p=2, m=2, p_int_first=0.1, p_int_others=0.01, p_repeat=0.05, p_complex=0.5)[source]¶
Create a discrete-time, random, stable, linear state space model.
- Parameters:
n (int, optional) – Number of states. Default is 2.
p (int, optional) – Number of control inputs. Default is 2.
m (int, optional) – Number of output measurements. If m=0, C becomes the identity matrix, so that y=x. Default is 2.
p_int_first (float, optional) – Probability of an integrator as the first pole. Default is 0.1.
p_int_others (float, optional) – Probability of other integrators beyond the first. Default is 0.01.
p_repeat (float, optional) – Probability of repeated roots. Default is 0.05.
p_complex (float, optional) – Probability of complex roots. Default is 0.5.
- Returns:
A tuple containing the state transition matrix (A), control matrix (B), and measurement matrix (C).
A (numpy.ndarray): State transition matrix of shape (n, n). B (numpy.ndarray): Control matrix of shape (n, p). C (numpy.ndarray): Measurement matrix of shape (m, n). If m = 0, C is the
identity matrix.
- Return type:
Tuple[numpy.ndarray, numpy.ndarray, numpy.ndarray]
- pykoopman.common.advance_linear_system(x0, u, n, A=None, B=None, C=None)[source]¶
Simulate the linear system dynamics for a given number of steps.
- Parameters:
x0 (numpy.ndarray) – Initial state vector of shape (n,).
u (numpy.ndarray) – Control input array of shape (p,) or (p, n-1). If 1-dimensional, it will be converted to a row vector.
n (int) – Number of steps to simulate.
A (numpy.ndarray, optional) – State transition matrix of shape (n, n). If not provided, it defaults to None.
B (numpy.ndarray, optional) – Control matrix of shape (n, p). If not provided, it defaults to None.
C (numpy.ndarray, optional) – Measurement matrix of shape (m, n). If not provided, it defaults to None.
- Returns:
A tuple containing the state trajectory (x) and the output trajectory (y).
x (numpy.ndarray): State trajectory of shape (n, len(x0)). y (numpy.ndarray): Output trajectory of shape (n, C.shape[0]).
- Return type:
Tuple[numpy.ndarray, numpy.ndarray]
- class pykoopman.common.torus_dynamics(n_states=128, sparsity=5, freq_max=15, noisemag=0.0)[source]¶
Bases:
object
Initializes a torus_dynamics object.
- Parameters:
n_states (int, optional) – Number of states. Default is 128.
sparsity (int, optional) – Degree of sparsity. Default is 5.
freq_max (int, optional) – Maximum frequency. Default is 15.
noisemag (float, optional) – Magnitude of noise. Default is 0.0.
- advance(n_samples, dt=1)[source]¶
Advances the continuous-time dynamics without control.
- Parameters:
n_samples (int) – Number of samples to advance.
dt (float, optional) – Time step. Default is 1.
- advance_discrete_time(n_samples, dt, u=None)[source]¶
Advances the discrete-time dynamics with or without control.
- Parameters:
n_samples (int) – Number of samples to advance.
dt (float) – Time step.
u (array-like, optional) – Control input. Default is None.
- set_control_matrix_physical(B)[source]¶
Sets the control matrix in physical space.
- Parameters:
B (array-like) – Control matrix in physical space.
- set_control_matrix_fourier(Bhat)[source]¶
Sets the control matrix in Fourier space.
- Parameters:
Bhat (array-like) – Control matrix in Fourier space.
- set_point_actuator(position=None)[source]¶
Sets a single point actuator.
- Parameters:
position (array-like, optional) – Position of the actuator. Default is None.
- viz_torus(ax, x)[source]¶
Visualizes the torus dynamics.
- Parameters:
ax – Axes object for plotting.
x (array-like) – Dynamics to be visualized.
- Returns:
Surface plot of the torus dynamics.
- Return type:
surface
- viz_all_modes(modes=None)[source]¶
Visualizes all modes.
- Parameters:
modes (array-like, optional) – Modes to be visualized. Default is None.
- Returns:
Figure object containing the visualizations.
- Return type:
fig
- property modes¶
Returns the modes of the dynamics.
- Returns:
Modes of the dynamics.
- Return type:
modes (array-like)
- property B_effective¶
Returns the effective control matrix.
- Returns:
Effective control matrix.
- Return type:
B_effective (array-like)
- pykoopman.common.lorenz(x, t, sigma=10, beta=2.6666666666666665, rho=28)[source]¶
Compute the derivative of the Lorenz system at a given state.
- Parameters:
x (list) – Current state of the Lorenz system [x, y, z].
t (float) – Current time.
sigma (float, optional) – Parameter sigma. Default is 10.
beta (float, optional) – Parameter beta. Default is 8/3.
rho (float, optional) – Parameter rho. Default is 28.
- Returns:
Derivative of the Lorenz system [dx/dt, dy/dt, dz/dt].
- Return type:
list
- pykoopman.common.vdp_osc(t, x, u)[source]¶
Compute the dynamics of the Van der Pol oscillator.
- Parameters:
t (float) – Time.
x (numpy.ndarray) – State vector of shape (2,).
u (float) – Control input.
- Returns:
Updated state vector of shape (2,).
- Return type:
numpy.ndarray
- pykoopman.common.rk4(t, x, u, _dt=0.01, func=<function vdp_osc>)[source]¶
Perform a 4th order Runge-Kutta integration.
- Parameters:
t (float) – Time.
x (numpy.ndarray) – State vector of shape (2,).
u (float) – Control input.
_dt (float, optional) – Time step. Defaults to 0.01.
func (function, optional) – Function defining the dynamics. Defaults to vdp_osc.
- Returns:
Updated state vector of shape (2,).
- Return type:
numpy.ndarray
- pykoopman.common.rev_dvdp(t, x, u=0, dt=0.1)[source]¶
Reverse dynamics of the Van der Pol oscillator.
- Parameters:
t (float) – Time.
x (numpy.ndarray) – Current state of the system [x1, x2].
u (float, optional) – Input. Default is 0.
dt (float, optional) – Time step. Default is 0.1.
- Returns:
Updated state of the system [x1’, x2’].
- Return type:
numpy.ndarray
- class pykoopman.common.Linear2Ddynamics[source]¶
Bases:
object
Initializes a Linear2Ddynamics object.
- linear_map(x)[source]¶
Applies the linear mapping to the input state.
- Parameters:
x (numpy.ndarray) – Input state.
- Returns:
Resulting mapped state.
- Return type:
numpy.ndarray
- class pykoopman.common.slow_manifold(mu=-0.05, la=-1.0, dt=0.01)[source]¶
Bases:
object
Represents the slow manifold class.
- Parameters:
mu (float, optional) – Parameter mu. Default is -0.05.
la (float, optional) – Parameter la. Default is -1.0.
dt (float, optional) – Time step size. Default is 0.01.
- Attributes:
mu (float) – Parameter mu.
la (float) – Parameter la.
b (float) – Value computed from mu and la.
dt (float) – Time step size.
n_states (int) – Number of states.
- sys(t, x, u)[source]¶
Computes the system dynamics.
- Parameters:
t (float) – Time.
x (array-like) – State.
u (array-like) – Control input.
- Returns:
Computed system dynamics.
- Return type:
array-like
- output(x)[source]¶
Computes the output based on the state.
- Parameters:
x (array-like) – State.
- Returns:
Computed output.
- Return type:
array-like
- simulate(x0, n_int)[source]¶
Simulates the system dynamics.
- Parameters:
x0 (array-like) – Initial state.
n_int (int) – Number of integration steps.
- Returns:
Simulated trajectory.
- Return type:
array-like
- collect_data_continuous(x0)[source]¶
Collects data from continuous-time dynamics.
- Parameters:
x0 (array-like) – Initial state.
- Returns:
Collected data (X, Y).
- Return type:
tuple
- collect_data_discrete(x0, n_int)[source]¶
Collects data from discrete-time dynamics.
- Parameters:
x0 (array-like) – Initial state.
n_int (int) – Number of integration steps.
- Returns:
Collected data (X, Y).
- Return type:
tuple
- class pykoopman.common.nlse(n, dt, L=6.283185307179586)[source]¶
Bases:
object
nonlinear schrodinger equation
iu_t + 0.5u_xx + u*|u|^2 = 0
periodic B.C. PDE is solved with Spectral methods using FFT
- collect_data_continuous(x0)[source]¶
collect training data pairs - continuous sense.
given x0, with shape (n_dim, n_traj), the function returns dx/dt with shape (n_dim, n_traj)
- class pykoopman.common.vbe(n, x, dt, nu=0.1, L=6.283185307179586)[source]¶
Bases:
object
1D viscous Burgers equation
u_t = -u*u_x +
u u_{xx}
periodic B.C. PDE is solved using spectral methods
- collect_data_continuous(x0)[source]¶
collect training data pairs - continuous sense.
given x0, with shape (n_dim, n_traj), the function returns dx/dt with shape (n_dim, n_traj)
- class pykoopman.common.cqgle(n, x, dt, tau=0.08, kappa=0, beta=0.66, nu=-0.1, sigma=-0.1, gamma=-0.1, L=6.283185307179586)[source]¶
Bases:
object
Cubic-quintic Ginzburg-Landau equation solver.
Solves the equation: i*u_t + (0.5 - i * tau) u_{xx} - i * kappa u_{xxxx} + (1-i * beta)|u|^2 u + (nu - i * sigma)|u|^4 u - i * gamma u = 0
Solves the periodic boundary conditions PDE using spectral methods.
- Attributes:
n_states (int) – Number of states.
x (numpy.ndarray) – x-coordinates.
dt (float) – Time step.
tau (float) – Parameter tau.
kappa (float) – Parameter kappa.
beta (float) – Parameter beta.
nu (float) – Parameter nu.
sigma (float) – Parameter sigma.
gamma (float) – Parameter gamma.
k (numpy.ndarray) – Wave numbers.
dk (float) – Wave number spacing.
- collect_data_continuous(x0)[source]¶
collect training data pairs - continuous sense.
given x0, with shape (n_dim, n_traj), the function returns dx/dt with shape (n_dim, n_traj)
- class pykoopman.common.ks(n, x, nu, dt, M=16)[source]¶
Bases:
object
Solving 1D KS equation
u_t = -u*u_x + u_{xx} + nu*u_{xxxx}
Periodic B.C. between 0 and 2*pi. This PDE is solved using spectral methods.