NystroemKernelRidge¶
- class kooplearn.kernel.NystroemKernelRidge(n_components=None, *, lag_time=1, reduced_rank=True, kernel='linear', gamma=None, degree=3, coef0=1, kernel_params=None, alpha=1e-06, eigen_solver='auto', tol=0, max_iter=None, n_centers=0.1, random_state=None, copy_X=True, n_jobs=None)[source]¶
Bases:
BaseEstimatorNyström-accelerated Kernel model minimizing the \(L^{2}\) loss.
Implements a model approximating the Koopman (deterministic systems) or Transfer (stochastic systems) operator by lifting the state with a infinite-dimensional nonlinear feature map associated to a kernel \(k\) and then minimizing the \(L^{2}\) loss in the embedded space as described in Meanti et al. [1].
Tip
The dynamical modes obtained by calling
kooplearn.kernel.NystroemKernelRidge.dynamical_modescorrespond to the Kernel Dynamical Mode Decomposition by O. Williams et al. [2].- Parameters:
n_components (
intorNone, defaultNone) – Number of components to retain. IfNone, all components are used.lag_time (
int, default1) – Time delay between the pairs of snapshots \((X_t, X_{t + \text{lag_time}})\) used to train the estimator.reduced_rank (
bool, defaultTrue) – Whether to use reduced-rank regression introduced in Kostic et al. [3]. IfFalse, initializes the classical principal component estimator.kernel (
{'linear', 'poly', 'rbf', 'sigmoid', 'cosine'}) – or callable, default=’linear’ Kernel function to use, or a callable that returns a Gram matrix.gamma (
floatorNone, defaultNone) – Kernel coefficient forrbf,poly, andsigmoidkernels. Ignored by other kernels. IfNone, it is set to1 / n_features.degree (
float, default3) – Degree for polynomial kernels. Ignored by other kernels.coef0 (
float, default1) – Independent term in polynomial and sigmoid kernels. Ignored by other kernels.kernel_params (
dictorNone, defaultNone) – Parameters (keyword arguments) passed to a callable kernel. Ignored by other kernels.alpha (
float, default1e-6) – Tikhonov (ridge) regularization coefficient.Noneis equivalent toalpha = 0, and internally calls specialized stable algorithms to deal with this specific case.eigen_solver (
{'auto', 'dense', 'arpack'}, default'auto') –Solver used to perform the internal SVD calcuations. If
n_componentsis much less than the number of training samples,arpackmay be more efficient than thedenseeigensolver.- auto :
the solver is selected automatically based on the number of samples and components: if the number of components to extract is less than 10 (strict) and the number of samples is more than 200 (strict), the
arpackmethod is enabled. Otherwise the exact full eigenvalue decomposition is computed and optionally truncated afterwards (densemethod).- dense :
run exact full eigenvalue decomposition calling the standard LAPACK solver via
scipy.linalg.eigh, and select the components by postprocessing.- arpack :
run SVD truncated to
n_componentscalling ARPACK solver usingscipy.sparse.linalg.eigsh. It requires strictly0 < n_components < n_samples.
tol (
float, default0) – Convergence tolerance for arpack. If 0, optimal value is chosen automatically by arpack.max_iter (
intorNone, defaultNone) – Maximum number of iterations for arpack. IfNone, optimal value is chosen automatically by arpack.n_centers (
intorfloat, default0.1) – Number of centers to select for the Nyström approximation. Ifn_centers < 1, selectsint(n_centers * n_samples)centers. Ifn_centers >= 1, selectsint(n_centers)centers.random_state (
int,RandomState instanceorNone, defaultNone) – Used wheneigen_solver == 'arpack'. Pass an int for reproducible results across multiple function calls.copy_X (
bool, defaultTrue) – IfTrue, input X is copied and stored by the model in theX_fit_attribute. If no further changes will be done to X, settingcopy_X=Falsesaves memory by storing a reference.n_jobs (
intorNone, defaultNone) – Number of parallel jobs to run.Nonemeans 1 unless in ajoblib.parallel_backendcontext.-1uses all available processors.
- Variables:
nystroem_centers_idxs (
ndarrayofint) – Indices of the selected centers.X_fit (
ndarrayofshape (n_samples,n_features)) – The data used to fit the model. Ifcopy_X=False, thenX_fit_is a reference to the original data.y_fit (
ndarrayofshape (n_samples,)orNone) – The observable used to fit the model. If no observable is provided during fitting, this attribute isNone.gamma (
float) – Effective kernel coefficient for RBF, polynomial, and sigmoid kernels. Whengammais explicitly provided, this is the same asgamma. WhengammaisNone, this is the inferred value.kernel_X (
ndarrayofshape (n_centers,n_centers)) – Kernel matrix evaluated at the initial states.kernel_Y (
ndarrayofshape (n_centers,n_centers)) – Kernel matrix evaluated at the evolved states.kernel_YX (
ndarrayofshape (n_centers,n_centers)) – Cross-kernel matrix between evolved and initial states.U (
ndarrayofshape (n_centers,rank)) – Left projection matrix of the operator approximation: \(k(\cdot, X) U V^\top k(\cdot, Y)\) (see Kostic et al. [3]).V (
ndarrayofshape (n_centers,rank)) – Right projection matrix of the operator approximation: \(k(\cdot, X) U V^\top k(\cdot, Y)\) (see Kostic et al. [3]).rank (
int) – Effective rank of the fitted estimator.
Examples
>>> from kooplearn.datasets import make_linear_system >>> from kooplearn.kernel import NystroemKernelRidge >>> import numpy as np >>> >>> # Generate a linear system >>> A = np.array([[0.9, 0.1], [-0.1, 0.9]]) >>> X0 = np.array([1.0, 0.0]) >>> data = make_linear_system(X0, A, n_steps=100, noise=0.1, random_state=42).to_numpy() >>> >>> # Fit the model >>> model = NystroemKernelRidge(n_components=2, kernel='linear', alpha=1e-3) >>> model = model.fit(data) >>> >>> # Predict the future state >>> pred = model.predict(data) >>> >>> # Get the eigenvalues of the Koopman operator >>> eigvals = model.eig()
Methods
- dynamical_modes(X, observable=False) DynamicalModes[source]¶
Compute the mode decomposition of arbitrary observables of the evolution operator at the states defined by
X. If \((\lambda_i, \xi_i, \psi_i)_{i = 1}^{r}\) are eigentriplets of the evolution operator, for any observable \(f\) the i-th mode of \(f\) at \(x\) is defined as: \(\lambda_i \langle \xi_i, f \rangle \psi_i(x)\). See for more details.- Parameters:
X (
ndarrayofshape (n_samples,n_features)) – States at which to evaluate the modes.observable (
bool, defaultFalse) – IfTrue, computes the modes of the observable rather than those of the system state.
- Returns:
Object containing the eigenvalues, modes, and their projections. See
kooplearn.structs.DynamicalModesfor details.- Return type:
DynamicalModes
- eig(eval_left_on=None, eval_right_on=None)[source]¶
Compute the eigendecomposition of the learned evolution operator.
This method returns the eigenvalues of the estimated evolution operator, and optionally evaluates the corresponding left and/or right eigenfunctions on user-provided data.
- Parameters:
eval_left_on (
ndarrayofshape (n_samples,n_features), optional) – Data points on which to evaluate the left eigenfunctions. IfNone, left eigenfunctions are not evaluated.eval_right_on (
ndarrayofshape (n_samples,n_features), optional) – Data points on which to evaluate the right eigenfunctions. IfNone, right eigenfunctions are not evaluated.
- Returns:
eigenvalues (
ndarrayofshape (n_components,)) – Eigenvalues of the estimated operator.left_eigenfunctions (
ndarrayofshape (n_samples,n_components), optional) – Values of the left eigenfunctions evaluated oneval_left_on. Returned only ifeval_left_onis provided.right_eigenfunctions (
ndarrayofshape (n_samples,n_components), optional) – Values of the right eigenfunctions evaluated oneval_right_on. Returned only ifeval_right_onis provided.
- fit(X, y=None)[source]¶
Fit the Nyström Kernel model.
Depending on the model parameters, this method estimates the evolution operator using a Nyström approximation of the kernel feature space using either a full-rank or reduced-rank regression approach.
- Parameters:
X (
ndarrayofshape (n_samples,n_features)) – Training trajectory data.y (
ndarrayofshape (n_samples,n_features_out), defaultNone) – Optional observable used for training. IfNone, the observable is assumed to be the state itself.
- Returns:
self – Fitted model instance.
- Return type:
object
- predict(X, n_steps=1, observable=False) ndarray[source]¶
Predict the system state or its expected value after
n_steps.Computes the predicted state — or, in the case of a stochastic system, the expected value \(\mathbb{E}[X_t \mid X_0 = X]\) — after
t = n_stepstime steps given the initial conditionsX. Ifobservable=True, returns the corresponding predicted observable instead of the state.- Parameters:
X (
ndarrayofshape (n_samples,n_features)) – Initial conditions used for prediction.n_steps (
int, default1) – Number of future time steps to predict. Only the final predicted state (at timet = n_steps) is returned.observable (
bool, defaultFalse) – IfTrue, returns the predicted observable at time \(t\) instead of the system state.
- Returns:
Predicted (expected) state or observable at time \(t = n_steps\).
- Return type:
ndarrayofshape (n_samples,n_features)
- risk(X=None)[source]¶
Compute the estimator risk.
- Parameters:
X (
ndarrayorNone, defaultNone) – Trajectory used to evaluate the risk. If None, evaluates on training data.- Returns:
Risk of the estimator.
- Return type:
float
- score(X=None, y=None, n_steps=1, observable=False, metric=<function r2_score>, **metric_kws) float[source]¶
Score the model predictions for timestep
n_steps.Computes the
metric(default issklearn.metrics.r2_score) evaluated between the model predictions at timestepn_stepsand the true system state (or observable, ifobservable=True) at the same timestep.- Parameters:
X – ndarray of shape (n_samples, n_features) or None, default=None Trajectory used to compute the score. If None, evaluates on training data.
y – ndarray of shape (n_samples, n_features) or None, default=None Optional observable used to compute the score.
n_steps – int, default=1 Number of future time steps on which to compute the score. Only the predictions at the final timestep (
t = n_steps) are compared to the true system state (or observable).observable – bool, default=False If
True, returns the predicted observable at time \(t\) instead of the system state.metric – callable (default=r2_score) The metric function used to score the model predictions.
metric_kws – dict Optional parameters to pass to the metric function.
- Returns:
- float
Metric function value for the model predictions at the next timestep.
- Return type:
score
Giacomo Meanti, Antoine Chatalic, Vladimir R Kostic, Pietro Novelli, Massimiliano Pontil, and Lorenzo Rosasco. Estimating koopman operators with sketching to provably learn large scale dynamical systems. In Thirty-seventh Conference on Neural Information Processing Systems. 2023. URL: https://openreview.net/forum?id=GItLpB1vhK.
Matthew O. Williams, Clarence W. Rowley, and Ioannis G. Kevrekidis. A kernel-based method for data-driven koopman spectral analysis. Journal of Computational Dynamics, 2(2):247–265, 2015. URL: http://dx.doi.org/10.3934/jcd.2015005, doi:10.3934/jcd.2015005.