QuadraticOperator#
- class QuadraticOperator(entries=None)[source]#
Quadratic state operator \(\Ophat_{\ell}(\qhat,\u) = \Hhat[\qhat\otimes\qhat]\) where \(\Hhat\in\RR^{r \times r^{2}}\).
Internally, the action of the operator is computed as the product of a \(r \times r(r+1)/2\) matrix and a compressed version of the Kronecker product \(\qhat \otimes \qhat\).
- Parameters
- entries(r, r^2) or (r, r(r+1)/2) or (r, r, r) ndarray or None
Operator entries \(\Hhat\).
Examples
>>> import numpy as np >>> H = opinf.operators.QuadraticOperator() >>> entries = np.random.random((10, 100)) # Operator entries. >>> H.set_entries(entries) >>> H.shape # Compressed shape. (10, 55) >>> q = np.random.random(10) # State vector. >>> out = H.apply(q) # Apply the operator to q. >>> np.allclose(out, entries @ np.kron(q, q)) True
Properties
entries
Discrete representation of the operator, the matrix \(\Ohat\).
shape
Shape of the operator entries array.
state_dimension
Dimension of the state \(\qhat\) that the operator acts on.
Methods
Apply the operator to the given state / input: \(\Ophat_{\ell}(\qhat,\u) = \Hhat[\qhat\otimes\qhat]\)
Calculate the compressed Kronecker product of a vector with itself.
Construct a mask for efficiently computing the compressed Kronecker product.
Given \(\Hhat\in\RR^{a\times r^2}\), construct the matrix \(\tilde{\H}\in\RR^{a \times r(r+1)/2}\) such that \(\Hhat[\qhat\otimes\qhat] = \tilde{\H}[\qhat\hat{\otimes}\qhat]\) for all \(\qhat\in\RR^{r}\) where \(\hat{\otimes}\) is the compressed Kronecker product (see
ckron()
).Return a copy of the operator.
Return the data matrix block corresponding to the operator, the Khatri--Rao product of the state with itself: \(\Qhat\odot\Qhat\) where \(\Qhat\) is
states
.Given \(\tilde{\H}\in\RR^{a \times r(r+1)/2}\), construct the matrix \(\Hhat\in\RR^{a\times r^2}\) such that \(\Hhat[\qhat\otimes\qhat] = \tilde{\H}[\qhat\hat{\otimes}\qhat]\) for all \(\qhat\in\RR^{r}\) where \(\hat{\otimes}\) is the compressed Kronecker product (see
ckron()
).Return the (Petrov-)Galerkin projection of the operator, \(\Hhat = \Wr\trp\H(\Vr\otimes\Vr)\).
Construct the state Jacobian of the operator: \(\ddqhat\Ophat_{\ell}(\qhat,\u) = \Hhat[(\I_r\otimes\qhat) + (\qhat\otimes\I_r)]\).
Load an operator from an HDF5 file.
Column dimension \(r(r+1)/2\) of the operator entries.
Save the operator to an HDF5 file.
Set the
entries
attribute.