pymor.models.iosys
¶
Module Contents¶
Classes¶
Base class for input-output systems. |
|
Base class for input-output systems with state space. |
|
Class for linear time-invariant systems. |
|
Class for systems represented by a transfer function. |
|
Class for linear second order systems. |
|
Class for linear delay systems. |
|
Class for linear stochastic systems. |
|
Class for bilinear systems. |
Functions¶
Return minimal sparse problem size for which to warn about converting to dense. |
|
Compute poles and residues. |
|
Create an |
- pymor.models.iosys.sparse_min_size(value=1000)[source]¶
Return minimal sparse problem size for which to warn about converting to dense.
- class pymor.models.iosys.InputOutputModel(dim_input, dim_output, cont_time=True, error_estimator=None, visualizer=None, name=None)[source]¶
Bases:
pymor.models.interface.Model
Base class for input-output systems.
- freq_resp(self, w, mu=None)[source]¶
Evaluate the transfer function on the imaginary axis.
Parameters
- w
A sequence of angular frequencies at which to compute the transfer function.
- mu
Parameter values
for which to evaluate the transfer function.
Returns
- tfw
Transfer function values at frequencies in
w
,NumPy array
of shape(len(w), self.dim_output, self.dim_input)
.
- bode(self, w, mu=None)[source]¶
Compute magnitudes and phases.
Parameters
- w
A sequence of angular frequencies at which to compute the transfer function.
- mu
Parameter values
for which to evaluate the transfer function.
Returns
- mag
Transfer function magnitudes at frequencies in
w
,NumPy array
of shape(len(w), self.dim_output, self.dim_input)
.- phase
Transfer function phases (in radians) at frequencies in
w
,NumPy array
of shape(len(w), self.dim_output, self.dim_input)
.
- bode_plot(self, w, mu=None, ax=None, Hz=False, dB=False, deg=True, **mpl_kwargs)[source]¶
Draw the Bode plot for all input-output pairs.
Parameters
- w
A sequence of angular frequencies at which to compute the transfer function.
- mu
Parameter
for which to evaluate the transfer function.- ax
Axis of shape (2 *
self.dim_output
,self.dim_input
) to which to plot. If not given,matplotlib.pyplot.gcf
is used to get the figure and create axis.- Hz
Should the frequency be in Hz on the plot.
- dB
Should the magnitude be in dB on the plot.
- deg
Should the phase be in degrees (otherwise in radians).
- mpl_kwargs
Keyword arguments used in the matplotlib plot function.
Returns
- artists
List of matplotlib artists added.
- mag_plot(self, w, mu=None, ax=None, ord=None, Hz=False, dB=False, **mpl_kwargs)[source]¶
Draw the magnitude plot.
Parameters
- w
A sequence of angular frequencies at which to compute the transfer function.
- mu
Parameter values
for which to evaluate the transfer function.- ax
Axis to which to plot. If not given,
matplotlib.pyplot.gca
is used.- ord
The order of the norm used to compute the magnitude (the default is the Frobenius norm).
- Hz
Should the frequency be in Hz on the plot.
- dB
Should the magnitude be in dB on the plot.
- mpl_kwargs
Keyword arguments used in the matplotlib plot function.
Returns
- out
List of matplotlib artists added.
- h2_norm(self, return_norm_only=True, **quad_kwargs)[source]¶
Compute the H2-norm using quadrature.
This method uses
scipy.integrate.quad
and makes no assumptions on the form of the transfer function.By default, the absolute error tolerance in
scipy.integrate.quad
is set to zero (see its optional argumentepsabs
). It can be changed by using theepsabs
keyword argument.Parameters
- return_norm_only
Whether to only return the approximate H2-norm.
- quad_kwargs
Keyword arguments passed to
scipy.integrate.quad
.
Returns
- norm
Computed H2-norm.
- norm_relerr
Relative error estimate (returned if
return_norm_only
isFalse
).- info
Quadrature info (returned if
return_norm_only
isFalse
andfull_output
isTrue
). Seescipy.integrate.quad
documentation for more details.
- h2_inner(self, lti)[source]¶
Compute H2 inner product with an
LTIModel
.Uses the inner product formula based on the pole-residue form (see, e.g., Lemma 1 in [ABG10]).
Parameters
- lti
LTIModel
consisting ofOperators
that can be converted toNumPy arrays
. The D operator is ignored.
Returns
- inner
H2 inner product.
- class pymor.models.iosys.InputStateOutputModel(dim_input, solution_space, dim_output, cont_time=True, error_estimator=None, visualizer=None, name=None)[source]¶
Bases:
InputOutputModel
Base class for input-output systems with state space.
- class pymor.models.iosys.LTIModel(A, B, C, D=None, E=None, cont_time=True, solver_options=None, error_estimator=None, visualizer=None, name=None)[source]¶
Bases:
InputStateOutputModel
Class for linear time-invariant systems.
This class describes input-state-output systems given by
\[\begin{split}E(\mu) \dot{x}(t, \mu) & = A(\mu) x(t, \mu) + B(\mu) u(t), \\ y(t, \mu) & = C(\mu) x(t, \mu) + D(\mu) u(t),\end{split}\]if continuous-time, or
\[\begin{split}E(\mu) x(k + 1, \mu) & = A(\mu) x(k, \mu) + B(\mu) u(k), \\ y(k, \mu) & = C(\mu) x(k, \mu) + D(\mu) u(k),\end{split}\]if discrete-time, where \(A\), \(B\), \(C\), \(D\), and \(E\) are linear operators.
Parameters
- A
The
Operator
A.- B
The
Operator
B.- C
The
Operator
C.- D
The
Operator
D orNone
(then D is assumed to be zero).- E
The
Operator
E orNone
(then E is assumed to be identity).- cont_time
True
if the system is continuous-time, otherwiseFalse
.- solver_options
The solver options to use to solve the Lyapunov equations.
- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)
method. Iferror_estimator
is notNone
, anestimate_error(U, mu)
method is added to the model which will callerror_estimator.estimate_error(U, mu, self)
.- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)
method. Ifvisualizer
is notNone
, avisualize(U, *args, **kwargs)
method is added to the model which forwards its arguments to the visualizer’svisualize
method.- name
Name of the system.
- classmethod from_matrices(cls, A, B, C, D=None, E=None, cont_time=True, state_id='STATE', solver_options=None, error_estimator=None, visualizer=None, name=None)[source]¶
Create
LTIModel
from matrices.Parameters
- A
The
NumPy array
orSciPy spmatrix
A.- B
The
NumPy array
orSciPy spmatrix
B.- C
The
NumPy array
orSciPy spmatrix
C.- D
The
NumPy array
orSciPy spmatrix
D orNone
(then D is assumed to be zero).- E
The
NumPy array
orSciPy spmatrix
E orNone
(then E is assumed to be identity).- cont_time
True
if the system is continuous-time, otherwiseFalse
.- state_id
Id of the state space.
- solver_options
The solver options to use to solve the Lyapunov equations.
- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)
method. Iferror_estimator
is notNone
, anestimate_error(U, mu)
method is added to the model which will callerror_estimator.estimate_error(U, mu, self)
.- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)
method. Ifvisualizer
is notNone
, avisualize(U, *args, **kwargs)
method is added to the model which forwards its arguments to the visualizer’svisualize
method.- name
Name of the system.
Returns
- lti
The
LTIModel
with operators A, B, C, D, and E.
- to_matrices(self)[source]¶
Return operators as matrices.
Returns
- A
The
NumPy array
orSciPy spmatrix
A.- B
The
NumPy array
orSciPy spmatrix
B.- C
The
NumPy array
orSciPy spmatrix
C.- D
The
NumPy array
orSciPy spmatrix
D orNone
(if D is aZeroOperator
).- E
The
NumPy array
orSciPy spmatrix
E orNone
(if E is anIdentityOperator
).
- classmethod from_files(cls, A_file, B_file, C_file, D_file=None, E_file=None, cont_time=True, state_id='STATE', solver_options=None, error_estimator=None, visualizer=None, name=None)[source]¶
Create
LTIModel
from matrices stored in separate files.Parameters
- A_file
The name of the file (with extension) containing A.
- B_file
The name of the file (with extension) containing B.
- C_file
The name of the file (with extension) containing C.
- D_file
None
or the name of the file (with extension) containing D.- E_file
None
or the name of the file (with extension) containing E.- cont_time
True
if the system is continuous-time, otherwiseFalse
.- state_id
Id of the state space.
- solver_options
The solver options to use to solve the Lyapunov equations.
- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)
method. Iferror_estimator
is notNone
, anestimate_error(U, mu)
method is added to the model which will callerror_estimator.estimate_error(U, mu, self)
.- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)
method. Ifvisualizer
is notNone
, avisualize(U, *args, **kwargs)
method is added to the model which forwards its arguments to the visualizer’svisualize
method.- name
Name of the system.
Returns
- lti
The
LTIModel
with operators A, B, C, D, and E.
- to_files(self, A_file, B_file, C_file, D_file=None, E_file=None)[source]¶
Write operators to files as matrices.
Parameters
- A_file
The name of the file (with extension) containing A.
- B_file
The name of the file (with extension) containing B.
- C_file
The name of the file (with extension) containing C.
- D_file
The name of the file (with extension) containing D or
None
if D is aZeroOperator
.- E_file
The name of the file (with extension) containing E or
None
if E is anIdentityOperator
.
- classmethod from_mat_file(cls, file_name, cont_time=True, state_id='STATE', solver_options=None, error_estimator=None, visualizer=None, name=None)[source]¶
Create
LTIModel
from matrices stored in a .mat file.Parameters
- file_name
The name of the .mat file (extension .mat does not need to be included) containing A, B, C, and optionally D and E.
- cont_time
True
if the system is continuous-time, otherwiseFalse
.- state_id
Id of the state space.
- solver_options
The solver options to use to solve the Lyapunov equations.
- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)
method. Iferror_estimator
is notNone
, anestimate_error(U, mu)
method is added to the model which will callerror_estimator.estimate_error(U, mu, self)
.- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)
method. Ifvisualizer
is notNone
, avisualize(U, *args, **kwargs)
method is added to the model which forwards its arguments to the visualizer’svisualize
method.- name
Name of the system.
Returns
- lti
The
LTIModel
with operators A, B, C, D, and E.
- to_mat_file(self, file_name)[source]¶
Save operators as matrices to .mat file.
Parameters
- file_name
The name of the .mat file (extension .mat does not need to be included).
- classmethod from_abcde_files(cls, files_basename, cont_time=True, state_id='STATE', solver_options=None, error_estimator=None, visualizer=None, name=None)[source]¶
Create
LTIModel
from matrices stored in .[ABCDE] files.Parameters
- files_basename
The basename of files containing A, B, C, and optionally D and E.
- cont_time
True
if the system is continuous-time, otherwiseFalse
.- state_id
Id of the state space.
- solver_options
The solver options to use to solve the Lyapunov equations.
- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)
method. Iferror_estimator
is notNone
, anestimate_error(U, mu)
method is added to the model which will callerror_estimator.estimate_error(U, mu, self)
.- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)
method. Ifvisualizer
is notNone
, avisualize(U, *args, **kwargs)
method is added to the model which forwards its arguments to the visualizer’svisualize
method.- name
Name of the system.
Returns
- lti
The
LTIModel
with operators A, B, C, D, and E.
- to_abcde_files(self, files_basename)[source]¶
Save operators as matrices to .[ABCDE] files in Matrix Market format.
Parameters
- files_basename
The basename of files containing the operators.
- poles(self, mu=None)[source]¶
Compute system poles.
Note
Assumes the systems is small enough to use a dense eigenvalue solver.
Parameters
- mu
Parameter values
for which to compute the systems poles.
Returns
One-dimensional
NumPy array
of system poles.
- eval_tf(self, s, mu=None)[source]¶
Evaluate the transfer function.
The transfer function at \(s\) is
\[C(\mu) (s E(\mu) - A(\mu))^{-1} B(\mu) + D(\mu).\]Note
Assumes that the number of inputs and outputs is much smaller than the order of the system.
Parameters
- s
Complex number.
- mu
Returns
- tfs
Transfer function evaluated at the complex number
s
,NumPy array
of shape(self.dim_output, self.dim_input)
.
- eval_dtf(self, s, mu=None)[source]¶
Evaluate the derivative of the transfer function.
The derivative of the transfer function at \(s\) is
\[-C(\mu) (s E(\mu) - A(\mu))^{-1} E(\mu) (s E(\mu) - A(\mu))^{-1} B(\mu).\]Note
Assumes that the number of inputs and outputs is much smaller than the order of the system.
Parameters
- s
Complex number.
- mu
Returns
- dtfs
Derivative of transfer function evaluated at the complex number
s
,NumPy array
of shape(self.dim_output, self.dim_input)
.
- gramian(self, typ, mu=None)[source]¶
Compute a Gramian.
Parameters
- typ
The type of the Gramian:
'c_lrcf'
: low-rank Cholesky factor of the controllability Gramian,'o_lrcf'
: low-rank Cholesky factor of the observability Gramian,'c_dense'
: dense controllability Gramian,'o_dense'
: dense observability Gramian.
Note
For
'c_lrcf'
and'o_lrcf'
types, the method assumes the system is asymptotically stable. For'c_dense'
and'o_dense'
types, the method assumes there are no two system poles which add to zero.- mu
Returns
If typ is
'c_lrcf'
or'o_lrcf'
, then the Gramian factor as aVectorArray
fromself.A.source
. If typ is'c_dense'
or'o_dense'
, then the Gramian as aNumPy array
.
- _hsv_U_V(self, mu=None)[source]¶
Compute Hankel singular values and vectors.
Note
Assumes the system is asymptotically stable.
Parameters
Returns
- hsv
One-dimensional
NumPy array
of singular values.- Uh
NumPy array
of left singular vectors.- Vh
NumPy array
of right singular vectors.
- hsv(self, mu=None)[source]¶
Hankel singular values.
Note
Assumes the system is asymptotically stable.
Parameters
Returns
- sv
One-dimensional
NumPy array
of singular values.
- h2_norm(self, mu=None)[source]¶
Compute the H2-norm of the
LTIModel
.Note
Assumes the system is asymptotically stable.
Parameters
Returns
- norm
H_2-norm.
- hinf_norm(self, mu=None, return_fpeak=False, ab13dd_equilibrate=False)[source]¶
Compute the H_infinity-norm of the
LTIModel
.Note
Assumes the system is asymptotically stable. Under this is assumption the H_infinity-norm is equal to the L_infinity-norm. Accordingly, this method calls
linf_norm
.Parameters
- mu
- return_fpeak
Whether to return the frequency at which the maximum is achieved.
- ab13dd_equilibrate
Whether
slycot.ab13dd
should use equilibration.
Returns
- norm
H_infinity-norm.
- fpeak
Frequency at which the maximum is achieved (if
return_fpeak
isTrue
).
- hankel_norm(self, mu=None)[source]¶
Compute the Hankel-norm of the
LTIModel
.Note
Assumes the system is asymptotically stable.
Parameters
Returns
- norm
Hankel-norm.
- l2_norm(self, ast_pole_data=None, mu=None)[source]¶
Compute the L2-norm of the
LTIModel
.The L2-norm of an
LTIModel
is defined via the integral\[\lVert H \rVert_{\mathcal{L}_2} = \left( \frac{1}{2 \pi} \int_{-\infty}^{\infty} \lVert H(\boldsymbol{\imath} \omega) \rVert_{\operatorname{F}}^2 \operatorname{d}\!\omega \right)^{\frac{1}{2}}.\]Parameters
- ast_pole_data
Can be:
dictionary of parameters for
eigs
,list of anti-stable eigenvalues (scalars),
tuple
(lev, ew, rev)
whereew
contains the anti-stable eigenvalues andlev
andrev
areVectorArrays
representing the eigenvectors.None
if anti-stable eigenvalues should be computed via dense methods.
- mu
Returns
- norm
L_2-norm.
- linf_norm(self, mu=None, return_fpeak=False, ab13dd_equilibrate=False)[source]¶
Compute the L_infinity-norm of the
LTIModel
.The L-infinity norm of an
LTIModel
is defined via\[\lVert H \rVert_{\mathcal{L}_\infty} = \sup_{\omega \in \mathbb{R}} \lVert H(\boldsymbol{\imath} \omega) \rVert_2.\]Parameters
- mu
- return_fpeak
Whether to return the frequency at which the maximum is achieved.
- ab13dd_equilibrate
Whether
slycot.ab13dd
should use equilibration.
Returns
- norm
L_infinity-norm.
- fpeak
Frequency at which the maximum is achieved (if
return_fpeak
isTrue
).
- get_ast_spectrum(self, ast_pole_data=None, mu=None)[source]¶
Compute anti-stable subset of the poles of the
LTIModel
.Parameters
- ast_pole_data
Can be:
dictionary of parameters for
eigs
,list of anti-stable eigenvalues (scalars),
tuple
(lev, ew, rev)
whereew
contains the sorted anti-stable eigenvalues andlev
andrev
areVectorArrays
representing the eigenvectors.None
if anti-stable eigenvalues should be computed via dense methods.
- mu
Returns
- lev
VectorArray
of left eigenvectors.- ew
One-dimensional
NumPy array
of anti-stable eigenvalues sorted from smallest to largest.- rev
VectorArray
of right eigenvectors.
- class pymor.models.iosys.TransferFunction(dim_input, dim_output, tf, dtf, parameters={}, cont_time=True, name=None)[source]¶
Bases:
InputOutputModel
Class for systems represented by a transfer function.
This class describes input-output systems given by a transfer function \(H(s, mu)\).
Parameters
- dim_input
The number of inputs.
- dim_output
The number of outputs.
- tf
The transfer function defined at least on the open right complex half-plane.
tf(s, mu)
is aNumPy array
of shape(p, m)
.- dtf
The complex derivative of
H
with respect tos
.- cont_time
True
if the system is continuous-time, otherwiseFalse
.- name
Name of the system.
- class pymor.models.iosys.SecondOrderModel(M, E, K, B, Cp, Cv=None, D=None, cont_time=True, solver_options=None, error_estimator=None, visualizer=None, name=None)[source]¶
Bases:
InputStateOutputModel
Class for linear second order systems.
This class describes input-output systems given by
\[\begin{split}M(\mu) \ddot{x}(t, \mu) + E(\mu) \dot{x}(t, \mu) + K(\mu) x(t, \mu) & = B(\mu) u(t), \\ y(t, \mu) & = C_p(\mu) x(t, \mu) + C_v(\mu) \dot{x}(t, \mu) + D(\mu) u(t),\end{split}\]if continuous-time, or
\[\begin{split}M(\mu) x(k + 2, \mu) + E(\mu) x(k + 1, \mu) + K(\mu) x(k, \mu) & = B(\mu) u(k), \\ y(k, \mu) & = C_p(\mu) x(k, \mu) + C_v(\mu) x(k + 1, \mu) + D(\mu) u(k),\end{split}\]if discrete-time, where \(M\), \(E\), \(K\), \(B\), \(C_p\), \(C_v\), and \(D\) are linear operators.
Parameters
- M
The
Operator
M.- E
The
Operator
E.- K
The
Operator
K.- B
The
Operator
B.- Cp
The
Operator
Cp.- Cv
The
Operator
Cv orNone
(then Cv is assumed to be zero).- D
The
Operator
D orNone
(then D is assumed to be zero).- cont_time
True
if the system is continuous-time, otherwiseFalse
.- solver_options
The solver options to use to solve the Lyapunov equations.
- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)
method. Iferror_estimator
is notNone
, anestimate_error(U, mu)
method is added to the model which will callerror_estimator.estimate_error(U, mu, self)
.- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)
method. Ifvisualizer
is notNone
, avisualize(U, *args, **kwargs)
method is added to the model which forwards its arguments to the visualizer’svisualize
method.- name
Name of the system.
- classmethod from_matrices(cls, M, E, K, B, Cp, Cv=None, D=None, cont_time=True, state_id='STATE', solver_options=None, error_estimator=None, visualizer=None, name=None)[source]¶
Create a second order system from matrices.
Parameters
- M
The
NumPy array
orSciPy spmatrix
M.- E
The
NumPy array
orSciPy spmatrix
E.- K
The
NumPy array
orSciPy spmatrix
K.- B
The
NumPy array
orSciPy spmatrix
B.- Cp
The
NumPy array
orSciPy spmatrix
Cp.- Cv
The
NumPy array
orSciPy spmatrix
Cv orNone
(then Cv is assumed to be zero).- D
The
NumPy array
orSciPy spmatrix
D orNone
(then D is assumed to be zero).- cont_time
True
if the system is continuous-time, otherwiseFalse
.- solver_options
The solver options to use to solve the Lyapunov equations.
- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)
method. Iferror_estimator
is notNone
, anestimate_error(U, mu)
method is added to the model which will callerror_estimator.estimate_error(U, mu, self)
.- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)
method. Ifvisualizer
is notNone
, avisualize(U, *args, **kwargs)
method is added to the model which forwards its arguments to the visualizer’svisualize
method.- name
Name of the system.
Returns
- lti
The SecondOrderModel with operators M, E, K, B, Cp, Cv, and D.
- to_matrices(self)[source]¶
Return operators as matrices.
Returns
- M
The
NumPy array
orSciPy spmatrix
M.- E
The
NumPy array
orSciPy spmatrix
E.- K
The
NumPy array
orSciPy spmatrix
K.- B
The
NumPy array
orSciPy spmatrix
B.- Cp
The
NumPy array
orSciPy spmatrix
Cp.- Cv
The
NumPy array
orSciPy spmatrix
Cv orNone
(if Cv is aZeroOperator
).- D
The
NumPy array
orSciPy spmatrix
D orNone
(if D is aZeroOperator
).
- classmethod from_files(cls, M_file, E_file, K_file, B_file, Cp_file, Cv_file=None, D_file=None, cont_time=True, state_id='STATE', solver_options=None, error_estimator=None, visualizer=None, name=None)[source]¶
Create
LTIModel
from matrices stored in separate files.Parameters
- M_file
The name of the file (with extension) containing A.
- E_file
The name of the file (with extension) containing E.
- K_file
The name of the file (with extension) containing K.
- B_file
The name of the file (with extension) containing B.
- Cp_file
The name of the file (with extension) containing Cp.
- Cv_file
None
or the name of the file (with extension) containing Cv.- D_file
None
or the name of the file (with extension) containing D.- cont_time
True
if the system is continuous-time, otherwiseFalse
.- state_id
Id of the state space.
- solver_options
The solver options to use to solve the Lyapunov equations.
- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)
method. Iferror_estimator
is notNone
, anestimate_error(U, mu)
method is added to the model which will callerror_estimator.estimate_error(U, mu, self)
.- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)
method. Ifvisualizer
is notNone
, avisualize(U, *args, **kwargs)
method is added to the model which forwards its arguments to the visualizer’svisualize
method.- name
Name of the system.
Returns
- some
The
SecondOrderModel
with operators M, E, K, B, Cp, Cv, and D.
- to_files(self, M_file, E_file, K_file, B_file, Cp_file, Cv_file=None, D_file=None)[source]¶
Write operators to files as matrices.
Parameters
- M_file
The name of the file (with extension) containing M.
- E_file
The name of the file (with extension) containing E.
- K_file
The name of the file (with extension) containing K.
- B_file
The name of the file (with extension) containing B.
- Cp_file
The name of the file (with extension) containing Cp.
- Cv_file
The name of the file (with extension) containing Cv or
None
if D is aZeroOperator
.- D_file
The name of the file (with extension) containing D or
None
if D is aZeroOperator
.
- to_lti(self)[source]¶
Return a first order representation.
The first order representation
\[\begin{split}\begin{bmatrix} I & 0 \\ 0 & M \end{bmatrix} \frac{\mathrm{d}}{\mathrm{d}t}\! \begin{bmatrix} x(t) \\ \dot{x}(t) \end{bmatrix} & = \begin{bmatrix} 0 & I \\ -K & -E \end{bmatrix} \begin{bmatrix} x(t) \\ \dot{x}(t) \end{bmatrix} + \begin{bmatrix} 0 \\ B \end{bmatrix} u(t), \\ y(t) & = \begin{bmatrix} C_p & C_v \end{bmatrix} \begin{bmatrix} x(t) \\ \dot{x}(t) \end{bmatrix} + D u(t)\end{split}\]is returned.
Returns
- lti
LTIModel
equivalent to the second-order model.
- __add__(self, other)[source]¶
Add a
SecondOrderModel
or anLTIModel
.
- __sub__(self, other)[source]¶
Subtract a
SecondOrderModel
or anLTIModel
.
- __neg__(self)[source]¶
Negate the
SecondOrderModel
.
- __mul__(self, other)[source]¶
Postmultiply by a
SecondOrderModel
or anLTIModel
.
- poles(self, mu=None)[source]¶
Compute system poles.
Note
Assumes the systems is small enough to use a dense eigenvalue solver.
Parameters
Returns
One-dimensional
NumPy array
of system poles.
- eval_tf(self, s, mu=None)[source]¶
Evaluate the transfer function.
The transfer function at \(s\) is
\[(C_p(\mu) + s C_v(\mu)) (s^2 M(\mu) + s E(\mu) + K(\mu))^{-1} B(\mu) + D(\mu).\]Note
Assumes that the number of inputs and outputs is much smaller than the order of the system.
Parameters
- s
Complex number.
- mu
Returns
- tfs
Transfer function evaluated at the complex number
s
,NumPy array
of shape(self.dim_output, self.dim_input)
.
- eval_dtf(self, s, mu=None)[source]¶
Evaluate the derivative of the transfer function.
\[s C_v(\mu) (s^2 M(\mu) + s E(\mu) + K(\mu))^{-1} B(\mu) - (C_p(\mu) + s C_v(\mu)) (s^2 M(\mu) + s E(\mu) + K(\mu))^{-1} (2 s M(\mu) + E(\mu)) (s^2 M(\mu) + s E(\mu) + K(\mu))^{-1} B(\mu).\]Note
Assumes that the number of inputs and outputs is much smaller than the order of the system.
Parameters
- s
Complex number.
- mu
Returns
- dtfs
Derivative of transfer function evaluated at the complex number
s
,NumPy array
of shape(self.dim_output, self.dim_input)
.
- gramian(self, typ, mu=None)[source]¶
Compute a second-order Gramian.
Parameters
- typ
The type of the Gramian:
'pc_lrcf'
: low-rank Cholesky factor of the position controllability Gramian,'vc_lrcf'
: low-rank Cholesky factor of the velocity controllability Gramian,'po_lrcf'
: low-rank Cholesky factor of the position observability Gramian,'vo_lrcf'
: low-rank Cholesky factor of the velocity observability Gramian,'pc_dense'
: dense position controllability Gramian,'vc_dense'
: dense velocity controllability Gramian,'po_dense'
: dense position observability Gramian,'vo_dense'
: dense velocity observability Gramian.
Note
For
'*_lrcf'
types, the method assumes the system is asymptotically stable. For'*_dense'
types, the method assumes there are no two system poles which add to zero.- mu
Returns
If typ is
'pc_lrcf'
,'vc_lrcf'
,'po_lrcf'
or'vo_lrcf'
, then the Gramian factor as aVectorArray
fromself.M.source
. If typ is'pc_dense'
,'vc_dense'
,'po_dense'
or'vo_dense'
, then the Gramian as aNumPy array
.
- psv(self, mu=None)[source]¶
Position singular values.
Note
Assumes the system is asymptotically stable.
Parameters
Returns
One-dimensional
NumPy array
of singular values.
- vsv(self, mu=None)[source]¶
Velocity singular values.
Note
Assumes the system is asymptotically stable.
Parameters
Returns
One-dimensional
NumPy array
of singular values.
- pvsv(self, mu=None)[source]¶
Position-velocity singular values.
Note
Assumes the system is asymptotically stable.
Parameters
Returns
One-dimensional
NumPy array
of singular values.
- vpsv(self, mu=None)[source]¶
Velocity-position singular values.
Note
Assumes the system is asymptotically stable.
Parameters
Returns
One-dimensional
NumPy array
of singular values.
- h2_norm(self, mu=None)[source]¶
Compute the H2-norm.
Note
Assumes the system is asymptotically stable.
Parameters
Returns
- norm
H_2-norm.
- hinf_norm(self, mu=None, return_fpeak=False, ab13dd_equilibrate=False)[source]¶
Compute the H_infinity-norm.
Note
Assumes the system is asymptotically stable.
Parameters
- mu
- return_fpeak
Should the frequency at which the maximum is achieved should be returned.
- ab13dd_equilibrate
Should
slycot.ab13dd
use equilibration.
Returns
- norm
H_infinity-norm.
- fpeak
Frequency at which the maximum is achieved (if
return_fpeak
isTrue
).
- class pymor.models.iosys.LinearDelayModel(A, Ad, tau, B, C, D=None, E=None, cont_time=True, error_estimator=None, visualizer=None, name=None)[source]¶
Bases:
InputStateOutputModel
Class for linear delay systems.
This class describes input-state-output systems given by
\[\begin{split}E x'(t) & = A x(t) + \sum_{i = 1}^q{A_i x(t - \tau_i)} + B u(t), \\ y(t) & = C x(t) + D u(t),\end{split}\]if continuous-time, or
\[\begin{split}E x(k + 1) & = A x(k) + \sum_{i = 1}^q{A_i x(k - \tau_i)} + B u(k), \\ y(k) & = C x(k) + D u(k),\end{split}\]if discrete-time, where \(E\), \(A\), \(A_i\), \(B\), \(C\), and \(D\) are linear operators.
Parameters
- A
The
Operator
A.- Ad
The tuple of
Operators
A_i.- tau
The tuple of delay times (positive floats or ints).
- B
The
Operator
B.- C
The
Operator
C.- D
The
Operator
D orNone
(then D is assumed to be zero).- E
The
Operator
E orNone
(then E is assumed to be identity).- cont_time
True
if the system is continuous-time, otherwiseFalse
.- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)
method. Iferror_estimator
is notNone
, anestimate_error(U, mu)
method is added to the model which will callerror_estimator.estimate_error(U, mu, self)
.- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)
method. Ifvisualizer
is notNone
, avisualize(U, *args, **kwargs)
method is added to the model which forwards its arguments to the visualizer’svisualize
method.- name
Name of the system.
- __add__(self, other)[source]¶
Add an
LTIModel
,SecondOrderModel
orLinearDelayModel
.
- __radd__(self, other)[source]¶
Add to an
LTIModel
or aSecondOrderModel
.
- __sub__(self, other)[source]¶
Subtract an
LTIModel
,SecondOrderModel
orLinearDelayModel
.
- __rsub__(self, other)[source]¶
Subtract from an
LTIModel
or aSecondOrderModel
.
- __neg__(self)[source]¶
Negate the
LinearDelayModel
.
- __mul__(self, other)[source]¶
Postmultiply an
LTIModel
,SecondOrderModel
orLinearDelayModel
.
- __rmul__(self, other)[source]¶
Premultiply by an
LTIModel
or aSecondOrderModel
.
- eval_tf(self, s, mu=None)[source]¶
Evaluate the transfer function.
The transfer function at \(s\) is
\[C \left(s E - A - \sum_{i = 1}^q{e^{-\tau_i s} A_i}\right)^{-1} B + D.\]Note
Assumes that the number of inputs and outputs is much smaller than the order of the system.
Parameters
- s
Complex number.
- mu
Returns
- tfs
Transfer function evaluated at the complex number
s
,NumPy array
of shape(self.dim_output, self.dim_input)
.
- eval_dtf(self, s, mu=None)[source]¶
Evaluate the derivative of the transfer function.
The derivative of the transfer function at \(s\) is
\[-C \left(s E - A - \sum_{i = 1}^q{e^{-\tau_i s} A_i}\right)^{-1} \left(E + \sum_{i = 1}^q{\tau_i e^{-\tau_i s} A_i}\right) \left(s E - A - \sum_{i = 1}^q{e^{-\tau_i s} A_i}\right)^{-1} B.\]Note
Assumes that the number of inputs and outputs is much smaller than the order of the system.
Parameters
- s
Complex number.
- mu
Returns
- dtfs
Derivative of transfer function evaluated at the complex number
s
,NumPy array
of shape(self.dim_output, self.dim_input)
.
- class pymor.models.iosys.LinearStochasticModel(A, As, B, C, D=None, E=None, cont_time=True, error_estimator=None, visualizer=None, name=None)[source]¶
Bases:
InputStateOutputModel
Class for linear stochastic systems.
This class describes input-state-output systems given by
\[\begin{split}E \mathrm{d}x(t) & = A x(t) \mathrm{d}t + \sum_{i = 1}^q{A_i x(t) \mathrm{d}\omega_i(t)} + B u(t) \mathrm{d}t, \\ y(t) & = C x(t) + D u(t),\end{split}\]if continuous-time, or
\[\begin{split}E x(k + 1) & = A x(k) + \sum_{i = 1}^q{A_i x(k) \omega_i(k)} + B u(k), \\ y(k) & = C x(k) + D u(t),\end{split}\]if discrete-time, where \(E\), \(A\), \(A_i\), \(B\), \(C\), and \(D\) are linear operators and \(\omega_i\) are stochastic processes.
Parameters
- A
The
Operator
A.- As
The tuple of
Operators
A_i.- B
The
Operator
B.- C
The
Operator
C.- D
The
Operator
D orNone
(then D is assumed to be zero).- E
The
Operator
E orNone
(then E is assumed to be identity).- cont_time
True
if the system is continuous-time, otherwiseFalse
.- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)
method. Iferror_estimator
is notNone
, anestimate_error(U, mu)
method is added to the model which will callerror_estimator.estimate_error(U, mu, self)
.- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)
method. Ifvisualizer
is notNone
, avisualize(U, *args, **kwargs)
method is added to the model which forwards its arguments to the visualizer’svisualize
method.- name
Name of the system.
- class pymor.models.iosys.BilinearModel(A, N, B, C, D, E=None, cont_time=True, error_estimator=None, visualizer=None, name=None)[source]¶
Bases:
InputStateOutputModel
Class for bilinear systems.
This class describes input-output systems given by
\[\begin{split}E x'(t) & = A x(t) + \sum_{i = 1}^m{N_i x(t) u_i(t)} + B u(t), \\ y(t) & = C x(t) + D u(t),\end{split}\]if continuous-time, or
\[\begin{split}E x(k + 1) & = A x(k) + \sum_{i = 1}^m{N_i x(k) u_i(k)} + B u(k), \\ y(k) & = C x(k) + D u(t),\end{split}\]if discrete-time, where \(E\), \(A\), \(N_i\), \(B\), \(C\), and \(D\) are linear operators and \(m\) is the number of inputs.
Parameters
- A
The
Operator
A.- N
The tuple of
Operators
N_i.- B
The
Operator
B.- C
The
Operator
C.- D
The
Operator
D orNone
(then D is assumed to be zero).- E
The
Operator
E orNone
(then E is assumed to be identity).- cont_time
True
if the system is continuous-time, otherwiseFalse
.- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)
method. Iferror_estimator
is notNone
, anestimate_error(U, mu)
method is added to the model which will callerror_estimator.estimate_error(U, mu, self)
.- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)
method. Ifvisualizer
is notNone
, avisualize(U, *args, **kwargs)
method is added to the model which forwards its arguments to the visualizer’svisualize
method.- name
Name of the system.
- pymor.models.iosys._lti_to_poles_b_c(lti)[source]¶
Compute poles and residues.
Parameters
- lti
LTIModel
consisting ofOperators
that can be converted toNumPy arrays
. The D operator is ignored.
Returns
- poles
1D
NumPy array
of poles.- b
NumPy array
of shape(lti.order, lti.dim_input)
.- c
NumPy array
of shape(lti.order, lti.dim_output)
.
- pymor.models.iosys._poles_b_c_to_lti(poles, b, c)[source]¶
Create an
LTIModel
from poles and residue rank-1 factors.Returns an
LTIModel
with real matrices such that its transfer function is\[\sum_{i = 1}^r \frac{c_i b_i^T}{s - \lambda_i}\]where \(\lambda_i, b_i, c_i\) are the poles and residue rank-1 factors.
Parameters
- poles
Sequence of poles.
- b
NumPy array
of shape(rom.order, rom.dim_input)
.- c
NumPy array
of shape(rom.order, rom.dim_output)
.
Returns