pykoopman package¶
Subpackages¶
- pykoopman.analytics package
- Module contents
BaseAnalyzer
ModesSelectionPAD21
PrunedKoopman
PrunedKoopman.fit
PrunedKoopman.predict
PrunedKoopman.psi
PrunedKoopman.phi
PrunedKoopman.ur
PrunedKoopman.A
PrunedKoopman.B
PrunedKoopman.C
PrunedKoopman.W
PrunedKoopman.lamda
PrunedKoopman.lamda_array
PrunedKoopman.continuous_lamda_array
PrunedKoopman.fit
PrunedKoopman.predict
PrunedKoopman.psi
PrunedKoopman.phi
PrunedKoopman.ur
PrunedKoopman.A
PrunedKoopman.B
PrunedKoopman.C
PrunedKoopman.W
PrunedKoopman.lamda
PrunedKoopman.lamda_array
PrunedKoopman.continuous_lamda_array
- Module contents
- pykoopman.common package
- Submodules
- pykoopman.common.cqgle module
cqgle
cqgle.sys
cqgle.simulate
cqgle.collect_data_continuous
cqgle.collect_one_step_data_discrete
cqgle.collect_one_trajectory_data
cqgle.visualize_data
cqgle.visualize_state_space
cqgle.sys
cqgle.simulate
cqgle.collect_data_continuous
cqgle.collect_one_step_data_discrete
cqgle.collect_one_trajectory_data
cqgle.visualize_data
cqgle.visualize_state_space
- pykoopman.common.examples module
drss
advance_linear_system
vdp_osc
rk4
square_wave
sine_wave
lorenz
rev_dvdp
Linear2Ddynamics
torus_dynamics
torus_dynamics.setup
torus_dynamics.advance
torus_dynamics.advance_discrete_time
torus_dynamics.set_control_matrix_physical
torus_dynamics.set_control_matrix_fourier
torus_dynamics.set_point_actuator
torus_dynamics.viz_setup
torus_dynamics.viz_torus
torus_dynamics.viz_all_modes
torus_dynamics.modes
torus_dynamics.B_effective
slow_manifold
slow_manifold.sys
slow_manifold.output
slow_manifold.simulate
slow_manifold.collect_data_continuous
slow_manifold.collect_data_discrete
slow_manifold.visualize_trajectories
slow_manifold.visualize_state_space
slow_manifold.sys
slow_manifold.output
slow_manifold.simulate
slow_manifold.collect_data_continuous
slow_manifold.collect_data_discrete
slow_manifold.visualize_trajectories
slow_manifold.visualize_state_space
forced_duffing
- pykoopman.common.ks module
- pykoopman.common.nlse module
- pykoopman.common.validation module
- pykoopman.common.vbe module
- Module contents
check_array
drop_nan_rows
validate_input
drss
advance_linear_system
torus_dynamics
torus_dynamics.setup
torus_dynamics.advance
torus_dynamics.advance_discrete_time
torus_dynamics.set_control_matrix_physical
torus_dynamics.set_control_matrix_fourier
torus_dynamics.set_point_actuator
torus_dynamics.viz_setup
torus_dynamics.viz_torus
torus_dynamics.viz_all_modes
torus_dynamics.modes
torus_dynamics.B_effective
lorenz
vdp_osc
rk4
rev_dvdp
Linear2Ddynamics
slow_manifold
slow_manifold.sys
slow_manifold.output
slow_manifold.simulate
slow_manifold.collect_data_continuous
slow_manifold.collect_data_discrete
slow_manifold.visualize_trajectories
slow_manifold.visualize_state_space
slow_manifold.sys
slow_manifold.output
slow_manifold.simulate
slow_manifold.collect_data_continuous
slow_manifold.collect_data_discrete
slow_manifold.visualize_trajectories
slow_manifold.visualize_state_space
nlse
vbe
cqgle
cqgle.sys
cqgle.simulate
cqgle.collect_data_continuous
cqgle.collect_one_step_data_discrete
cqgle.collect_one_trajectory_data
cqgle.visualize_data
cqgle.visualize_state_space
cqgle.sys
cqgle.simulate
cqgle.collect_data_continuous
cqgle.collect_one_step_data_discrete
cqgle.collect_one_trajectory_data
cqgle.visualize_data
cqgle.visualize_state_space
ks
- pykoopman.differentiation package
- pykoopman.observables package
- pykoopman.regression package
Submodules¶
pykoopman.koopman module¶
module for discrete time Koopman class
- class pykoopman.koopman.Koopman(observables=None, regressor=None, quiet=False)[source]¶
Bases:
BaseEstimator
Constructor for the Koopman class.
- Parameters:
observables – observables object, optional (default:
pykoopman.observables.Identity
) Map(s) to apply to raw measurement data before estimating the Koopman operator. Must extendpykoopman.observables.BaseObservables
. The default option,pykoopman.observables.Identity
, leaves the input untouched.regressor – regressor object, optional (default:
DMD
) The regressor used to learn the Koopman operator from the observables.regressor
can either extend thepykoopman.regression.BaseRegressor
, orpydmd.DMDBase
. In the latter case, the pydmd object must have both afit
and apredict
method.quiet – boolean, optional (default: False) Whether or not warnings should be silenced during fitting.
- fit(x, y=None, u=None, dt=1)[source]¶
Fit the Koopman model by learning an approximate Koopman operator.
- Parameters:
x – numpy.ndarray, shape (n_samples, n_features) Measurement data to be fit. Each row should correspond to an example and each column a feature. If only x is provided, it is assumed that examples are equi-spaced in time (i.e., a uniform timestep is assumed).
y – numpy.ndarray, shape (n_samples, n_features), optional (default: None) Target measurement data to be fit, i.e., it is assumed y = fun(x). Each row should correspond to an example and each column a feature. The samples in x and y are generally not required to be consecutive and equi-spaced.
u – numpy.ndarray, shape (n_samples, n_control_features), optional (default: None) Control/actuation/external parameter data. Each row should correspond to one sample and each column a control variable or feature. The control variable may be the amplitude of an actuator or an external, time-varying parameter. It is assumed that samples in u occur at the time instances of the corresponding samples in x, e.g., x(t+1) = fun(x(t), u(t)).
dt – float, optional (default: 1) Time step between samples
- Returns:
returns a fit
Koopman
instance- Return type:
self
- predict(x, u=None)[source]¶
Predict the state one timestep in the future.
- Parameters:
x – numpy.ndarray, shape (n_samples, n_input_features) Current state.
u – numpy.ndarray, shape (n_samples, n_control_features), optional (default None) Time series of external actuation/control.
- Returns:
- numpy.ndarray, shape (n_samples, n_input_features)
Predicted state one timestep in the future.
- Return type:
y
- simulate(x0, u=None, n_steps=1)[source]¶
Simulate an initial state forward in time with the learned Koopman model.
- Parameters:
x0 – numpy.ndarray, shape (n_input_features,) or (n_consumed_samples + 1, n_input_features) Initial state from which to simulate. If using TimeDelay observables,
x0
should contain enough examples to compute all required time delays, i.e.,n_consumed_samples + 1
.u – numpy.ndarray, shape (n_samples, n_control_features), optional (default None) Time series of external actuation/control.
n_steps – int, optional (default 1) Number of forward steps to be simulated.
- Returns:
- numpy.ndarray, shape (n_steps, n_input_features)
Simulated states. Note that
y[0, :]
is one timestep ahead ofx0
.
- Return type:
y
- get_feature_names(input_features=None)[source]¶
Get the names of the individual features constituting the observables.
- Parameters:
input_features – list of string, length n_input_features, optional (default None) String names for input features, if available. By default, the names “x0”, “x1”, …, “xn_input_features” are used.
- Returns:
- list of string, length n_output_features
Output feature names.
- Return type:
output_feature_names
- phi(x_col)[source]¶
Compute the feature matrix phi(x) given
x_col
.- Parameters:
x_col – numpy.ndarray, shape (n_features, n_samples) State vectors to be evaluated for phi.
- Returns:
- numpy.ndarray, shape (n_samples, self.n_output_features_)
Value of phi evaluated at input
x_col
.
- Return type:
phi
- psi(x_col)[source]¶
Compute the Koopman psi(x) given
x_col
.- Parameters:
x_col – numpy.ndarray, shape (n_features, n_samples) State vectors to be evaluated for psi.
- Returns:
- numpy.ndarray, shape (n_samples, self.n_output_features_)
Value of psi evaluated at input
x_col
.
- Return type:
eigen_phi
- property A¶
Returns the state transition matrix
A
.The state transition matrix A satisfies y’ = Ay or y’ = Ay + Bu, respectively, where y = g(x) and y is a low-rank representation.
- property B¶
Returns the control matrix
B
.The control matrix (or vector) B satisfies y’ = Ay + Bu. y is the reduced system state.
- property C¶
Returns the measurement matrix (or vector) C.
The measurement matrix C satisfies x = C * phi_r.
- property W¶
Returns the Koopman modes.
- property lamda¶
Returns the discrete-time Koopman lambda obtained from spectral decomposition.
- property lamda_array¶
Returns the discrete-time Koopman lambda as an array.
- property continuous_lamda_array¶
Returns the continuous-time Koopman lambda as an array.
- property ur¶
Returns the projection matrix Ur.
- validity_check(t, x)[source]¶
Perform a validity check of eigenfunctions.
The validity check tests the linearity of eigenfunctions phi(x(t)) == phi(x(0)) * exp(lambda*t).
- Parameters:
t – numpy.ndarray, shape (n_samples,) Time vector.
x – numpy.ndarray, shape (n_samples, n_input_features) State vectors to be checked.
- Returns:
- list
Sorted indices of eigenfunctions based on linearity error.
- linearity_error: list
Linearity error for each eigenfunction.
- Return type:
efun_index
- score(x, y=None, cast_as_real=True, metric=<function r2_score>, **metric_kws)[source]¶
Score the model predictions for the next timestep.
- Parameters:
x – numpy.ndarray, shape (n_samples, n_input_features) State measurements.
y – numpy.ndarray, shape (n_samples, n_input_features), optional (default None). State measurements one timestep in the future.
cast_as_real – bool, optional (default True) Whether to take the real part of predictions when computing the score.
metric – callable, optional (default r2_score) The metric function used to score the model predictions.
metric_kws – dict, optional Optional parameters to pass to the metric function.
- Returns:
- float
Metric function value for the model predictions at the next timestep.
- Return type:
score
pykoopman.koopman_continuous module¶
module for continuous time Koopman class
- class pykoopman.koopman_continuous.KoopmanContinuous(observables=None, differentiator=Derivative(k=1, kind='finite_difference'), regressor=None)[source]¶
Bases:
Koopman
Continuous-time Koopman class.
- Parameters:
observables – Observables object, optional (default: pykoopman.observables.Identity) Map(s) to apply to raw measurement data before estimating the Koopman operator. Must extend pykoopman.observables.BaseObservables. The default option, pykoopman.observables.Identity, leaves the input untouched.
differentiator – Callable, optional (default: centered difference) Function used to compute numerical derivatives. The function must have the call signature differentiator(x, t), where x is a 2D numpy ndarray of shape (n_samples, n_features) and t is a 1D numpy ndarray of shape (n_samples,).
regressor – Regressor object, optional (default: DMD) The regressor used to learn the Koopman operator from the observables. regressor can either extend pykoopman.regression.BaseRegressor, or the pydmd.DMDBase class. In the latter case, the pydmd object must have both a fit and a predict method.
- predict(x, dt=0, u=None)[source]¶
Predict using continuous-time Koopman model.
- Parameters:
x – numpy.ndarray State measurements. Each row should correspond to the system state at some point in time.
dt – float, optional (default: 0) Time step between measurements. If specified, the prediction is made for the given time step in the future.
u – numpy.ndarray, optional (default: None) Control input/actuation data. Each row should correspond to one sample and each column a control variable or feature.
- Returns:
- numpy.ndarray
Predicted state using the continuous-time Koopman model. Each row corresponds to the predicted state for the corresponding row in x.
- Return type:
output
- simulate(x, t=0, u=None)[source]¶
Simulate continuous-time Koopman model.
- Parameters:
x – numpy.ndarray Initial state from which to simulate. Each row corresponds to the system state at some point in time.
t – float, optional (default: 0) Time at which to simulate the system. If specified, the simulation is performed for the given time.
u – numpy.ndarray, optional (default: None) Control input/actuation data. Each row should correspond to one sample and each column a control variable or feature.
- Returns:
- numpy.ndarray
Simulated states of the system. Each row corresponds to the simulated state at a specific time point.
- Return type:
output
Module contents¶
- class pykoopman.Koopman(observables=None, regressor=None, quiet=False)[source]¶
Bases:
BaseEstimator
Constructor for the Koopman class.
- Parameters:
observables – observables object, optional (default:
pykoopman.observables.Identity
) Map(s) to apply to raw measurement data before estimating the Koopman operator. Must extendpykoopman.observables.BaseObservables
. The default option,pykoopman.observables.Identity
, leaves the input untouched.regressor – regressor object, optional (default:
DMD
) The regressor used to learn the Koopman operator from the observables.regressor
can either extend thepykoopman.regression.BaseRegressor
, orpydmd.DMDBase
. In the latter case, the pydmd object must have both afit
and apredict
method.quiet – boolean, optional (default: False) Whether or not warnings should be silenced during fitting.
- fit(x, y=None, u=None, dt=1)[source]¶
Fit the Koopman model by learning an approximate Koopman operator.
- Parameters:
x – numpy.ndarray, shape (n_samples, n_features) Measurement data to be fit. Each row should correspond to an example and each column a feature. If only x is provided, it is assumed that examples are equi-spaced in time (i.e., a uniform timestep is assumed).
y – numpy.ndarray, shape (n_samples, n_features), optional (default: None) Target measurement data to be fit, i.e., it is assumed y = fun(x). Each row should correspond to an example and each column a feature. The samples in x and y are generally not required to be consecutive and equi-spaced.
u – numpy.ndarray, shape (n_samples, n_control_features), optional (default: None) Control/actuation/external parameter data. Each row should correspond to one sample and each column a control variable or feature. The control variable may be the amplitude of an actuator or an external, time-varying parameter. It is assumed that samples in u occur at the time instances of the corresponding samples in x, e.g., x(t+1) = fun(x(t), u(t)).
dt – float, optional (default: 1) Time step between samples
- Returns:
returns a fit
Koopman
instance- Return type:
self
- predict(x, u=None)[source]¶
Predict the state one timestep in the future.
- Parameters:
x – numpy.ndarray, shape (n_samples, n_input_features) Current state.
u – numpy.ndarray, shape (n_samples, n_control_features), optional (default None) Time series of external actuation/control.
- Returns:
- numpy.ndarray, shape (n_samples, n_input_features)
Predicted state one timestep in the future.
- Return type:
y
- simulate(x0, u=None, n_steps=1)[source]¶
Simulate an initial state forward in time with the learned Koopman model.
- Parameters:
x0 – numpy.ndarray, shape (n_input_features,) or (n_consumed_samples + 1, n_input_features) Initial state from which to simulate. If using TimeDelay observables,
x0
should contain enough examples to compute all required time delays, i.e.,n_consumed_samples + 1
.u – numpy.ndarray, shape (n_samples, n_control_features), optional (default None) Time series of external actuation/control.
n_steps – int, optional (default 1) Number of forward steps to be simulated.
- Returns:
- numpy.ndarray, shape (n_steps, n_input_features)
Simulated states. Note that
y[0, :]
is one timestep ahead ofx0
.
- Return type:
y
- get_feature_names(input_features=None)[source]¶
Get the names of the individual features constituting the observables.
- Parameters:
input_features – list of string, length n_input_features, optional (default None) String names for input features, if available. By default, the names “x0”, “x1”, …, “xn_input_features” are used.
- Returns:
- list of string, length n_output_features
Output feature names.
- Return type:
output_feature_names
- phi(x_col)[source]¶
Compute the feature matrix phi(x) given
x_col
.- Parameters:
x_col – numpy.ndarray, shape (n_features, n_samples) State vectors to be evaluated for phi.
- Returns:
- numpy.ndarray, shape (n_samples, self.n_output_features_)
Value of phi evaluated at input
x_col
.
- Return type:
phi
- psi(x_col)[source]¶
Compute the Koopman psi(x) given
x_col
.- Parameters:
x_col – numpy.ndarray, shape (n_features, n_samples) State vectors to be evaluated for psi.
- Returns:
- numpy.ndarray, shape (n_samples, self.n_output_features_)
Value of psi evaluated at input
x_col
.
- Return type:
eigen_phi
- property A¶
Returns the state transition matrix
A
.The state transition matrix A satisfies y’ = Ay or y’ = Ay + Bu, respectively, where y = g(x) and y is a low-rank representation.
- property B¶
Returns the control matrix
B
.The control matrix (or vector) B satisfies y’ = Ay + Bu. y is the reduced system state.
- property C¶
Returns the measurement matrix (or vector) C.
The measurement matrix C satisfies x = C * phi_r.
- property W¶
Returns the Koopman modes.
- property lamda¶
Returns the discrete-time Koopman lambda obtained from spectral decomposition.
- property lamda_array¶
Returns the discrete-time Koopman lambda as an array.
- property continuous_lamda_array¶
Returns the continuous-time Koopman lambda as an array.
- property ur¶
Returns the projection matrix Ur.
- validity_check(t, x)[source]¶
Perform a validity check of eigenfunctions.
The validity check tests the linearity of eigenfunctions phi(x(t)) == phi(x(0)) * exp(lambda*t).
- Parameters:
t – numpy.ndarray, shape (n_samples,) Time vector.
x – numpy.ndarray, shape (n_samples, n_input_features) State vectors to be checked.
- Returns:
- list
Sorted indices of eigenfunctions based on linearity error.
- linearity_error: list
Linearity error for each eigenfunction.
- Return type:
efun_index
- score(x, y=None, cast_as_real=True, metric=<function r2_score>, **metric_kws)[source]¶
Score the model predictions for the next timestep.
- Parameters:
x – numpy.ndarray, shape (n_samples, n_input_features) State measurements.
y – numpy.ndarray, shape (n_samples, n_input_features), optional (default None). State measurements one timestep in the future.
cast_as_real – bool, optional (default True) Whether to take the real part of predictions when computing the score.
metric – callable, optional (default r2_score) The metric function used to score the model predictions.
metric_kws – dict, optional Optional parameters to pass to the metric function.
- Returns:
- float
Metric function value for the model predictions at the next timestep.
- Return type:
score
- class pykoopman.KoopmanContinuous(observables=None, differentiator=Derivative(k=1, kind='finite_difference'), regressor=None)[source]¶
Bases:
Koopman
Continuous-time Koopman class.
- Parameters:
observables – Observables object, optional (default: pykoopman.observables.Identity) Map(s) to apply to raw measurement data before estimating the Koopman operator. Must extend pykoopman.observables.BaseObservables. The default option, pykoopman.observables.Identity, leaves the input untouched.
differentiator – Callable, optional (default: centered difference) Function used to compute numerical derivatives. The function must have the call signature differentiator(x, t), where x is a 2D numpy ndarray of shape (n_samples, n_features) and t is a 1D numpy ndarray of shape (n_samples,).
regressor – Regressor object, optional (default: DMD) The regressor used to learn the Koopman operator from the observables. regressor can either extend pykoopman.regression.BaseRegressor, or the pydmd.DMDBase class. In the latter case, the pydmd object must have both a fit and a predict method.
- predict(x, dt=0, u=None)[source]¶
Predict using continuous-time Koopman model.
- Parameters:
x – numpy.ndarray State measurements. Each row should correspond to the system state at some point in time.
dt – float, optional (default: 0) Time step between measurements. If specified, the prediction is made for the given time step in the future.
u – numpy.ndarray, optional (default: None) Control input/actuation data. Each row should correspond to one sample and each column a control variable or feature.
- Returns:
- numpy.ndarray
Predicted state using the continuous-time Koopman model. Each row corresponds to the predicted state for the corresponding row in x.
- Return type:
output
- simulate(x, t=0, u=None)[source]¶
Simulate continuous-time Koopman model.
- Parameters:
x – numpy.ndarray Initial state from which to simulate. Each row corresponds to the system state at some point in time.
t – float, optional (default: 0) Time at which to simulate the system. If specified, the simulation is performed for the given time.
u – numpy.ndarray, optional (default: None) Control input/actuation data. Each row should correspond to one sample and each column a control variable or feature.
- Returns:
- numpy.ndarray
Simulated states of the system. Each row corresponds to the simulated state at a specific time point.
- Return type:
output