StateInputOperator#

class StateInputOperator(entries=None)[source]#

Linear state / input interaction operator \(\Ophat_{\ell}(\qhat,\u) = \Nhat[\u\otimes\qhat]\) where \(\Nhat \in \RR^{r \times rm}\).

Parameters:
entries(r, rm) ndarray or None

Operator matrix \(\Nhat\).

Examples

>>> import numpy as np
>>> N = opinf.operators.StateInputOperator()
>>> entries = np.random.random((10, 3))
>>> N.set_entries(entries)
>>> N.shape
(10, 3)
>>> q = np.random.random(10)                # State vector.
>>> u = np.random.random(3)                 # Input vector.
>>> out = N.apply(q, u)                     # Apply the operator to (q,u).
>>> np.allclose(out, entries @ np.kron(u, q))
True
Properties:
entries#

Operator matrix \(\Nhat\).

input_dimension#

Dimension \(m\) of the input \(\u\) that the operator acts on.

shape#

Shape \((r, rm)\) of the operator matrix \(\Nhat\).

state_dimension#

Dimension \(r\) of the state \(\qhat\) that the operator acts on.

Methods:

apply

Apply the operator to the given state / input: \(\Ophat_{\ell}(\qhat,\u) = \Nhat[\u\otimes\qhat]\).

copy

Return a copy of the operator.

datablock

Return the data matrix block corresponding to the operator, the Khatri--Rao product \(\U\odot\Qhat\) where \(\Qhat\) is states and \(\U\) is inputs.

galerkin

Return the Galerkin projection of the operator, \(\Nhat = (\Wr\trp\Vr)^{-1}\Wr\trp\N[\I_{m}\otimes\Vr]\).

jacobian

Construct the state Jacobian of the operator: \(\ddqhat\Ophat_{\ell}(\qhat,\u) = \sum_{i=1}^{m}u_{i}\Nhat_{i}\) where \(\Nhat=[~\Nhat_{1}~~\cdots~~\Nhat_{m}~]\) and each \(\Nhat_i\in\RR^{r\times r},~i=1,\ldots,m\).

load

Load an operator from an HDF5 file.

operator_dimension

Column dimension \(rm\) of the operator matrix \(\Nhat\).

save

Save the operator to an HDF5 file.

set_entries

Set the operator matrix \(\Nhat\).

verify

Verify consistency between dimension properties and required methods.