pymor.models.iosys

Module Contents

Classes

LTIModel

Class for linear time-invariant systems.

SecondOrderModel

Class for linear second order systems.

LinearDelayModel

Class for linear delay systems.

LinearStochasticModel

Class for linear stochastic systems.

BilinearModel

Class for bilinear systems.

Functions

sparse_min_size

Return minimal sparse problem size for which to warn about converting to dense.

_lti_to_poles_b_c

Compute poles and residues.

_poles_b_c_to_lti

Create an LTIModel from poles and residue rank-1 factors.

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.LTIModel(A, B, C, D=None, E=None, cont_time=True, solver_options=None, error_estimator=None, visualizer=None, name=None)[source]

Bases: pymor.models.interface.Model

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.

All methods related to the transfer function (e.g., frequency response calculation and Bode plots) are attached to the transfer_function attribute.

Parameters

A

The Operator A.

B

The Operator B.

C

The Operator C.

D

The Operator D or None (then D is assumed to be zero).

E

The Operator E or None (then E is assumed to be identity).

cont_time

True if the system is continuous-time, otherwise False.

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. If error_estimator is not None, an estimate_error(U, mu) method is added to the model which will call error_estimator.estimate_error(U, mu, self).

visualizer

A visualizer for the problem. This can be any object with a visualize(U, model, ...) method. If visualizer is not None, a visualize(U, *args, **kwargs) method is added to the model which forwards its arguments to the visualizer’s visualize method.

name

Name of the system.

order[source]

The order of the system.

dim_input[source]

The number of inputs.

dim_output[source]

The number of outputs.

A[source]

The Operator A.

B[source]

The Operator B.

C[source]

The Operator C.

D[source]

The Operator D.

E[source]

The Operator E.

transfer_function[source]

The transfer function.

__str__(self)[source]

Return str(self).

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 or SciPy spmatrix A.

B

The NumPy array or SciPy spmatrix B.

C

The NumPy array or SciPy spmatrix C.

D

The NumPy array or SciPy spmatrix D or None (then D is assumed to be zero).

E

The NumPy array or SciPy spmatrix E or None (then E is assumed to be identity).

cont_time

True if the system is continuous-time, otherwise False.

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. If error_estimator is not None, an estimate_error(U, mu) method is added to the model which will call error_estimator.estimate_error(U, mu, self).

visualizer

A visualizer for the problem. This can be any object with a visualize(U, model, ...) method. If visualizer is not None, a visualize(U, *args, **kwargs) method is added to the model which forwards its arguments to the visualizer’s visualize 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 or SciPy spmatrix A.

B

The NumPy array or SciPy spmatrix B.

C

The NumPy array or SciPy spmatrix C.

D

The NumPy array or SciPy spmatrix D or None (if D is a ZeroOperator).

E

The NumPy array or SciPy spmatrix E or None (if E is an IdentityOperator).

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, otherwise False.

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. If error_estimator is not None, an estimate_error(U, mu) method is added to the model which will call error_estimator.estimate_error(U, mu, self).

visualizer

A visualizer for the problem. This can be any object with a visualize(U, model, ...) method. If visualizer is not None, a visualize(U, *args, **kwargs) method is added to the model which forwards its arguments to the visualizer’s visualize 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 a ZeroOperator.

E_file

The name of the file (with extension) containing E or None if E is an IdentityOperator.

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, otherwise False.

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. If error_estimator is not None, an estimate_error(U, mu) method is added to the model which will call error_estimator.estimate_error(U, mu, self).

visualizer

A visualizer for the problem. This can be any object with a visualize(U, model, ...) method. If visualizer is not None, a visualize(U, *args, **kwargs) method is added to the model which forwards its arguments to the visualizer’s visualize 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, otherwise False.

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. If error_estimator is not None, an estimate_error(U, mu) method is added to the model which will call error_estimator.estimate_error(U, mu, self).

visualizer

A visualizer for the problem. This can be any object with a visualize(U, model, ...) method. If visualizer is not None, a visualize(U, *args, **kwargs) method is added to the model which forwards its arguments to the visualizer’s visualize 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.

__add__(self, other)[source]

Add an LTIModel.

__sub__(self, other)[source]

Subtract an LTIModel.

__neg__(self)[source]

Negate the LTIModel.

__mul__(self, other)[source]

Postmultiply by an LTIModel.

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.

This function is an alias for transfer_function.eval_tf. See eval_tf for a detailed documentation.

eval_dtf(self, s, mu=None)[source]

Evaluate the derivative of the transfer function.

This function is an alias for transfer_function.eval_dtf. See eval_dtf for a detailed documentation.

freq_resp(self, w, mu=None)[source]

Evaluate the transfer function on the imaginary axis.

This function is an alias for transfer_function.freq_resp. See freq_resp for a detailed documentation.

bode(self, w, mu=None)[source]

Compute magnitudes and phases.

This function is an alias for transfer_function.bode. See bode for a detailed documentation.

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.

This function is an alias for transfer_function.bode_plot. See bode_plot for a detailed documentation.

mag_plot(self, w, mu=None, ax=None, ord=None, Hz=False, dB=False, **mpl_kwargs)[source]

Draw the magnitude plot.

This function is an alias for transfer_function.mag_plot. See mag_plot for a detailed documentation.

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

Parameter values.

Returns

If typ is 'c_lrcf' or 'o_lrcf', then the Gramian factor as a VectorArray from self.A.source. If typ is 'c_dense' or 'o_dense', then the Gramian as a NumPy array.

_hsv_U_V(self, mu=None)[source]

Compute Hankel singular values and vectors.

Note

Assumes the system is asymptotically stable.

Parameters

mu

Parameter values.

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

mu

Parameter values.

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

mu

Parameter values.

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

Parameter values.

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 is True).

hankel_norm(self, mu=None)[source]

Compute the Hankel-norm of the LTIModel.

Note

Assumes the system is asymptotically stable.

Parameters

mu

Parameter values.

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) where ew contains the anti-stable eigenvalues and lev and rev are VectorArrays representing the eigenvectors.

  • None if anti-stable eigenvalues should be computed via dense methods.

mu

Parameter.

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

Parameter.

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 is True).

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) where ew contains the sorted anti-stable eigenvalues and lev and rev are VectorArrays representing the eigenvectors.

  • None if anti-stable eigenvalues should be computed via dense methods.

mu

Parameter.

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.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: pymor.models.interface.Model

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.

All methods related to the transfer function (e.g., frequency response calculation and Bode plots) are attached to the transfer_function attribute.

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 or None (then Cv is assumed to be zero).

D

The Operator D or None (then D is assumed to be zero).

cont_time

True if the system is continuous-time, otherwise False.

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. If error_estimator is not None, an estimate_error(U, mu) method is added to the model which will call error_estimator.estimate_error(U, mu, self).

visualizer

A visualizer for the problem. This can be any object with a visualize(U, model, ...) method. If visualizer is not None, a visualize(U, *args, **kwargs) method is added to the model which forwards its arguments to the visualizer’s visualize method.

name

Name of the system.

order[source]

The order of the system (equal to M.source.dim).

dim_input[source]

The number of inputs.

dim_output[source]

The number of outputs.

M[source]

The Operator M.

E[source]

The Operator E.

K[source]

The Operator K.

B[source]

The Operator B.

Cp[source]

The Operator Cp.

Cv[source]

The Operator Cv.

D[source]

The Operator D.

transfer_function[source]

The transfer function.

__str__(self)[source]

Return str(self).

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 or SciPy spmatrix M.

E

The NumPy array or SciPy spmatrix E.

K

The NumPy array or SciPy spmatrix K.

B

The NumPy array or SciPy spmatrix B.

Cp

The NumPy array or SciPy spmatrix Cp.

Cv

The NumPy array or SciPy spmatrix Cv or None (then Cv is assumed to be zero).

D

The NumPy array or SciPy spmatrix D or None (then D is assumed to be zero).

cont_time

True if the system is continuous-time, otherwise False.

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. If error_estimator is not None, an estimate_error(U, mu) method is added to the model which will call error_estimator.estimate_error(U, mu, self).

visualizer

A visualizer for the problem. This can be any object with a visualize(U, model, ...) method. If visualizer is not None, a visualize(U, *args, **kwargs) method is added to the model which forwards its arguments to the visualizer’s visualize 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 or SciPy spmatrix M.

E

The NumPy array or SciPy spmatrix E.

K

The NumPy array or SciPy spmatrix K.

B

The NumPy array or SciPy spmatrix B.

Cp

The NumPy array or SciPy spmatrix Cp.

Cv

The NumPy array or SciPy spmatrix Cv or None (if Cv is a ZeroOperator).

D

The NumPy array or SciPy spmatrix D or None (if D is a ZeroOperator).

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, otherwise False.

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. If error_estimator is not None, an estimate_error(U, mu) method is added to the model which will call error_estimator.estimate_error(U, mu, self).

visualizer

A visualizer for the problem. This can be any object with a visualize(U, model, ...) method. If visualizer is not None, a visualize(U, *args, **kwargs) method is added to the model which forwards its arguments to the visualizer’s visualize 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 a ZeroOperator.

D_file

The name of the file (with extension) containing D or None if D is a ZeroOperator.

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 an LTIModel.

__radd__(self, other)[source]

Add to an LTIModel.

__sub__(self, other)[source]

Subtract a SecondOrderModel or an LTIModel.

__rsub__(self, other)[source]

Subtract from an LTIModel.

__neg__(self)[source]

Negate the SecondOrderModel.

__mul__(self, other)[source]

Postmultiply by a SecondOrderModel or an LTIModel.

__rmul__(self, other)[source]

Premultiply by an LTIModel.

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.

Returns

One-dimensional NumPy array of system poles.

eval_tf(self, s, mu=None)[source]

Evaluate the transfer function.

This function is an alias for transfer_function.eval_tf. See eval_tf for a detailed documentation.

eval_dtf(self, s, mu=None)[source]

Evaluate the derivative of the transfer function.

This function is an alias for transfer_function.eval_dtf. See eval_dtf for a detailed documentation.

freq_resp(self, w, mu=None)[source]

Evaluate the transfer function on the imaginary axis.

This function is an alias for transfer_function.freq_resp. See freq_resp for a detailed documentation.

bode(self, w, mu=None)[source]

Compute magnitudes and phases.

This function is an alias for transfer_function.bode. See bode for a detailed documentation.

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.

This function is an alias for transfer_function.bode_plot. See bode_plot for a detailed documentation.

mag_plot(self, w, mu=None, ax=None, ord=None, Hz=False, dB=False, **mpl_kwargs)[source]

Draw the magnitude plot.

This function is an alias for transfer_function.mag_plot. See mag_plot for a detailed documentation.

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

Parameter values.

Returns

If typ is 'pc_lrcf', 'vc_lrcf', 'po_lrcf' or 'vo_lrcf', then the Gramian factor as a VectorArray from self.M.source. If typ is 'pc_dense', 'vc_dense', 'po_dense' or 'vo_dense', then the Gramian as a NumPy array.

psv(self, mu=None)[source]

Position singular values.

Note

Assumes the system is asymptotically stable.

Parameters

mu

Parameter values.

Returns

One-dimensional NumPy array of singular values.

vsv(self, mu=None)[source]

Velocity singular values.

Note

Assumes the system is asymptotically stable.

Parameters

mu

Parameter values.

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

mu

Parameter values.

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

mu

Parameter values.

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

mu

Parameter values.

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

Parameter values.

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 is True).

hankel_norm(self, mu=None)[source]

Compute the Hankel-norm.

Note

Assumes the system is asymptotically stable.

Parameters

mu

Parameter values.

Returns

norm

Hankel-norm.

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: pymor.models.interface.Model

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.

All methods related to the transfer function (e.g., frequency response calculation and Bode plots) are attached to the transfer_function attribute.

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 or None (then D is assumed to be zero).

E

The Operator E or None (then E is assumed to be identity).

cont_time

True if the system is continuous-time, otherwise False.

error_estimator

An error estimator for the problem. This can be any object with an estimate_error(U, mu, model) method. If error_estimator is not None, an estimate_error(U, mu) method is added to the model which will call error_estimator.estimate_error(U, mu, self).

visualizer

A visualizer for the problem. This can be any object with a visualize(U, model, ...) method. If visualizer is not None, a visualize(U, *args, **kwargs) method is added to the model which forwards its arguments to the visualizer’s visualize method.

name

Name of the system.

order[source]

The order of the system (equal to A.source.dim).

dim_input[source]

The number of inputs.

dim_output[source]

The number of outputs.

q[source]

The number of delay terms.

tau[source]

The tuple of delay times.

A[source]

The Operator A.

Ad[source]

The tuple of Operators A_i.

B[source]

The Operator B.

C[source]

The Operator C.

D[source]

The Operator D.

E[source]

The Operator E.

transfer_function[source]

The transfer function.

__str__(self)[source]

Return str(self).

__add__(self, other)[source]

Add an LTIModel, SecondOrderModel or LinearDelayModel.

__radd__(self, other)[source]

Add to an LTIModel or a SecondOrderModel.

__sub__(self, other)[source]

Subtract an LTIModel, SecondOrderModel or LinearDelayModel.

__rsub__(self, other)[source]

Subtract from an LTIModel or a SecondOrderModel.

__neg__(self)[source]

Negate the LinearDelayModel.

__mul__(self, other)[source]

Postmultiply an LTIModel, SecondOrderModel or LinearDelayModel.

__rmul__(self, other)[source]

Premultiply by an LTIModel or a SecondOrderModel.

eval_tf(self, s, mu=None)[source]

Evaluate the transfer function.

This function is an alias for transfer_function.eval_tf. See eval_tf for a detailed documentation.

eval_dtf(self, s, mu=None)[source]

Evaluate the derivative of the transfer function.

This function is an alias for transfer_function.eval_dtf. See eval_dtf for a detailed documentation.

freq_resp(self, w, mu=None)[source]

Evaluate the transfer function on the imaginary axis.

This function is an alias for transfer_function.freq_resp. See freq_resp for a detailed documentation.

bode(self, w, mu=None)[source]

Compute magnitudes and phases.

This function is an alias for transfer_function.bode. See bode for a detailed documentation.

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.

This function is an alias for transfer_function.bode_plot. See bode_plot for a detailed documentation.

mag_plot(self, w, mu=None, ax=None, ord=None, Hz=False, dB=False, **mpl_kwargs)[source]

Draw the magnitude plot.

This function is an alias for transfer_function.mag_plot. See mag_plot for a detailed documentation.

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: pymor.models.interface.Model

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 or None (then D is assumed to be zero).

E

The Operator E or None (then E is assumed to be identity).

cont_time

True if the system is continuous-time, otherwise False.

error_estimator

An error estimator for the problem. This can be any object with an estimate_error(U, mu, model) method. If error_estimator is not None, an estimate_error(U, mu) method is added to the model which will call error_estimator.estimate_error(U, mu, self).

visualizer

A visualizer for the problem. This can be any object with a visualize(U, model, ...) method. If visualizer is not None, a visualize(U, *args, **kwargs) method is added to the model which forwards its arguments to the visualizer’s visualize method.

name

Name of the system.

order[source]

The order of the system (equal to A.source.dim).

dim_input[source]

The number of inputs.

dim_output[source]

The number of outputs.

q[source]

The number of stochastic processes.

A[source]

The Operator A.

As[source]

The tuple of Operators A_i.

B[source]

The Operator B.

C[source]

The Operator C.

D[source]

The Operator D.

E[source]

The Operator E.

__str__(self)[source]

Return str(self).

class pymor.models.iosys.BilinearModel(A, N, B, C, D, E=None, cont_time=True, error_estimator=None, visualizer=None, name=None)[source]

Bases: pymor.models.interface.Model

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 or None (then D is assumed to be zero).

E

The Operator E or None (then E is assumed to be identity).

cont_time

True if the system is continuous-time, otherwise False.

error_estimator

An error estimator for the problem. This can be any object with an estimate_error(U, mu, model) method. If error_estimator is not None, an estimate_error(U, mu) method is added to the model which will call error_estimator.estimate_error(U, mu, self).

visualizer

A visualizer for the problem. This can be any object with a visualize(U, model, ...) method. If visualizer is not None, a visualize(U, *args, **kwargs) method is added to the model which forwards its arguments to the visualizer’s visualize method.

name

Name of the system.

order[source]

The order of the system (equal to A.source.dim).

dim_input[source]

The number of inputs.

dim_output[source]

The number of outputs.

A[source]

The Operator A.

N[source]

The tuple of Operators N_i.

B[source]

The Operator B.

C[source]

The Operator C.

D[source]

The Operator D.

E[source]

The Operator E.

__str__(self)[source]

Return str(self).

pymor.models.iosys._lti_to_poles_b_c(lti)[source]

Compute poles and residues.

Parameters

lti

LTIModel consisting of Operators that can be converted to NumPy 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

LTIModel.