pymor.models package

Submodules

basic module


class pymor.models.basic.InstationaryModel(*args, **kwargs)[source]

Bases: pymor.models.interface.Model

Generic class for models of instationary problems.

This class describes instationary problems given by the equations:

M * ∂_t u(t, μ) + L(u(μ), t, μ) = F(t, μ)
                        u(0, μ) = u_0(μ)

for t in [0,T], where L is a (possibly non-linear) time-dependent Operator, F is a time-dependent vector-like Operator, and u_0 the initial data. The mass Operator M is assumed to be linear.

Parameters

T

The final time T.

initial_data

The initial data u_0. Either a VectorArray of length 1 or (for the Parameter-dependent case) a vector-like Operator (i.e. a linear Operator with source.dim == 1) which applied to NumpyVectorArray(np.array([1])) will yield the initial data for given parameter values.

operator

The Operator L.

rhs

The right-hand side F.

mass

The mass Operator M. If None, the identity is assumed.

time_stepper

The time-stepper to be used by solve.

num_values

The number of returned vectors of the solution trajectory. If None, each intermediate vector that is calculated is returned.

output_functional

Operator mapping a given solution to the model output. In many applications, this will be a Functional, i.e. an Operator mapping to scalars. This is not required, however.

products

A dict of product Operators defined on the discrete space the problem is posed on. For each product with key 'x' a corresponding attribute x_product, as well as a norm method x_norm is added to the model.

estimator

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

visualizer

A visualizer for the problem. This can be any object with a visualize(U, m, ...) 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 model.

__str__()[source]

Return str(self).

to_lti()[source]

Convert model to LTIModel.

This method interprets the given model as an LTIModel in the following way:

- self.operator        -> A
self.rhs               -> B
self.output_functional -> C
None                   -> D
self.mass              -> E

class pymor.models.basic.StationaryModel(*args, **kwargs)[source]

Bases: pymor.models.interface.Model

Generic class for models of stationary problems.

This class describes discrete problems given by the equation:

L(u(μ), μ) = F(μ)

with a vector-like right-hand side F and a (possibly non-linear) operator L.

Note that even when solving a variational formulation where F is a functional and not a vector, F has to be specified as a vector-like Operator (mapping scalars to vectors). This ensures that in the complex case both L and F are anti-linear in the test variable.

Parameters

operator

The Operator L.

rhs

The vector F. Either a VectorArray of length 1 or a vector-like Operator.

output_functional

Operator mapping a given solution to the model output. In many applications, this will be a Functional, i.e. an Operator mapping to scalars. This is not required, however.

products

A dict of inner product Operators defined on the discrete space the problem is posed on. For each product with key 'x' a corresponding attribute x_product, as well as a norm method x_norm is added to the model.

estimator

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

visualizer

A visualizer for the problem. This can be any object with a visualize(U, m, ...) 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 model.

__str__()[source]

Return str(self).

interface module


class pymor.models.interface.Model(*args, **kwargs)[source]

Bases: pymor.core.cache.CacheableObject, pymor.parameters.base.ParametricObject

Interface for model objects.

A model object defines a discrete problem via its class and the Operators it contains. Furthermore, models can be solved for given parameter values resulting in a solution VectorArray.

solution_space

VectorSpace of the solution VectorArrays returned by solve.

output_space

VectorSpace of the model output VectorArrays returned by output (typically NumpyVectorSpace(k) where k is a small).

linear

True if the model describes a linear problem.

products

Dict of inner product operators associated with the model.

estimate(U, mu=None)[source]

Estimate the model error for a given solution.

The model error could be the error w.r.t. the analytical solution of the given problem or the model reduction error w.r.t. a corresponding high-dimensional Model.

Parameters

U

The solution obtained by solve.

mu

Parameter values for which U has been obtained.

Returns

The estimated error.

output(mu=None, **kwargs)[source]

Return the model output for given parameter values mu.

Parameters

mu

Parameter values for which to compute the output.

Returns

The computed model output as a VectorArray from output_space.

solve(mu=None, return_output=False, **kwargs)[source]

Solve the discrete problem for the parameter values mu.

The result will be cached in case caching has been activated for the given model.

Parameters

mu

Parameter values for which to solve.

return_output

If True, the model output for the given parameter values mu is returned as a VectorArray from output_space.

Returns

The solution VectorArray. When return_output is True, the output VectorArray is returned as second value.

visualize(U, **kwargs)[source]

Visualize a solution VectorArray U.

Parameters

U

The VectorArray from solution_space that shall be visualized.

kwargs

See docstring of self.visualizer.visualize.

iosys module


class pymor.models.iosys.BilinearModel(*args, **kwargs)[source]

Bases: pymor.models.iosys.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 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.

estimator

An error estimator for the problem. This can be any object with an estimate(U, mu, model) method. If estimator is not None, an estimate(U, mu) method is added to the model which will call estimator.estimate(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

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

input_dim

The number of inputs.

output_dim

The number of outputs.

A

The Operator A.

N

The tuple of Operators N_i.

B

The Operator B.

C

The Operator C.

D

The Operator D.

E

The Operator E.

__str__()[source]

Return str(self).


class pymor.models.iosys.InputOutputModel(*args, **kwargs)[source]

Bases: pymor.models.interface.Model

Base class for input-output systems.

eval_dtf(s, mu=None)[source]

Evaluate the derivative of the transfer function.

eval_tf(s, mu=None)[source]

Evaluate the transfer function.

freq_resp(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.output_dim, self.input_dim).

mag_plot(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.


class pymor.models.iosys.InputStateOutputModel(*args, **kwargs)[source]

Bases: pymor.models.iosys.InputOutputModel

Base class for input-output systems with state space.


class pymor.models.iosys.LTIModel(*args, **kwargs)[source]

Bases: pymor.models.iosys.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 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.

estimator

An error estimator for the problem. This can be any object with an estimate(U, mu, model) method. If estimator is not None, an estimate(U, mu) method is added to the model which will call estimator.estimate(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

The order of the system.

input_dim

The number of inputs.

output_dim

The number of outputs.

A

The Operator A.

B

The Operator B.

C

The Operator C.

D

The Operator D.

E

The Operator E.

__add__(other)[source]

Add an LTIModel.

__mul__(other)[source]

Postmultiply by an LTIModel.

__neg__()[source]

Negate the LTIModel.

__str__()[source]

Return str(self).

__sub__(other)[source]

Subtract an LTIModel.

eval_dtf(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 either the number of inputs or the number of outputs is much smaller than the order of the system.

Parameters

s

Complex number.

mu

Parameter values.

Returns

dtfs

Derivative of transfer function evaluated at the complex number s, NumPy array of shape (self.output_dim, self.input_dim).

eval_tf(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 either the number of inputs or the number of outputs is much smaller than the order of the system.

Parameters

s

Complex number.

mu

Parameter values.

Returns

tfs

Transfer function evaluated at the complex number s, NumPy array of shape (self.output_dim, self.input_dim).

classmethod from_abcde_files(files_basename, cont_time=True, state_id='STATE', solver_options=None, estimator=None, visualizer=None, name=None)[source]

Create LTIModel from matrices stored in a .[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.

estimator

An error estimator for the problem. This can be any object with an estimate(U, mu, model) method. If estimator is not None, an estimate(U, mu) method is added to the model which will call estimator.estimate(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.

classmethod from_files(A_file, B_file, C_file, D_file=None, E_file=None, cont_time=True, state_id='STATE', solver_options=None, 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.

estimator

An error estimator for the problem. This can be any object with an estimate(U, mu, model) method. If estimator is not None, an estimate(U, mu) method is added to the model which will call estimator.estimate(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.

classmethod from_mat_file(file_name, cont_time=True, state_id='STATE', solver_options=None, 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.

estimator

An error estimator for the problem. This can be any object with an estimate(U, mu, model) method. If estimator is not None, an estimate(U, mu) method is added to the model which will call estimator.estimate(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.

classmethod from_matrices(A, B, C, D=None, E=None, cont_time=True, state_id='STATE', solver_options=None, 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.

estimator

An error estimator for the problem. This can be any object with an estimate(U, mu, model) method. If estimator is not None, an estimate(U, mu) method is added to the model which will call estimator.estimate(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.

gramian(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.

h2_norm(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.

hankel_norm(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.

hinf_norm(mu=None, return_fpeak=False, ab13dd_equilibrate=False)[source]

Compute the H_infinity-norm of the LTIModel.

Note

Assumes the system is asymptotically stable.

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

hsv(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.

poles(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.


class pymor.models.iosys.LinearDelayModel(*args, **kwargs)[source]

Bases: pymor.models.iosys.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 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.

estimator

An error estimator for the problem. This can be any object with an estimate(U, mu, model) method. If estimator is not None, an estimate(U, mu) method is added to the model which will call estimator.estimate(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

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

input_dim

The number of inputs.

output_dim

The number of outputs.

q

The number of delay terms.

tau

The tuple of delay times.

A

The Operator A.

Ad

The tuple of Operators A_i.

B

The Operator B.

C

The Operator C.

D

The Operator D.

E

The Operator E.

__add__(other)[source]

Add an LTIModel, SecondOrderModel or LinearDelayModel.

__mul__(other)[source]

Postmultiply by a SecondOrderModel or an LTIModel.

__neg__()[source]

Negate the LinearDelayModel.

__radd__(other)[source]

Add to an LTIModel or a SecondOrderModel.

__rmul__(other)[source]

Premultiply by an LTIModel or a SecondOrderModel.

__rsub__(other)[source]

Subtract from an LTIModel or a SecondOrderModel.

__str__()[source]

Return str(self).

__sub__(other)[source]

Subtract an LTIModel, SecondOrderModel or LinearDelayModel.

eval_dtf(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 either the number of inputs or the number of outputs is much smaller than the order of the system.

Parameters

s

Complex number.

mu

Parameter values.

Returns

dtfs

Derivative of transfer function evaluated at the complex number s, NumPy array of shape (self.output_dim, self.input_dim).

eval_tf(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 either the number of inputs or the number of outputs is much smaller than the order of the system.

Parameters

s

Complex number.

mu

Parameter values.

Returns

tfs

Transfer function evaluated at the complex number s, NumPy array of shape (self.output_dim, self.input_dim).


class pymor.models.iosys.LinearStochasticModel(*args, **kwargs)[source]

Bases: pymor.models.iosys.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 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.

estimator

An error estimator for the problem. This can be any object with an estimate(U, mu, model) method. If estimator is not None, an estimate(U, mu) method is added to the model which will call estimator.estimate(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

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

input_dim

The number of inputs.

output_dim

The number of outputs.

q

The number of stochastic processes.

A

The Operator A.

As

The tuple of Operators A_i.

B

The Operator B.

C

The Operator C.

D

The Operator D.

E

The Operator E.

__str__()[source]

Return str(self).


class pymor.models.iosys.SecondOrderModel(*args, **kwargs)[source]

Bases: pymor.models.iosys.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 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.

estimator

An error estimator for the problem. This can be any object with an estimate(U, mu, model) method. If estimator is not None, an estimate(U, mu) method is added to the model which will call estimator.estimate(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

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

input_dim

The number of inputs.

output_dim

The number of outputs.

M

The Operator M.

E

The Operator E.

K

The Operator K.

B

The Operator B.

Cp

The Operator Cp.

Cv

The Operator Cv.

D

The Operator D.

__add__(other)[source]

Add a SecondOrderModel or an LTIModel.

__mul__(other)[source]

Postmultiply by a SecondOrderModel or an LTIModel.

__neg__()[source]

Negate the SecondOrderModel.

__radd__(other)[source]

Add to an LTIModel.

__rmul__(other)[source]

Premultiply by an LTIModel.

__rsub__(other)[source]

Subtract from an LTIModel.

__str__()[source]

Return str(self).

__sub__(other)[source]

Subtract a SecondOrderModel or an LTIModel.

eval_dtf(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 either the number of inputs or the number of outputs is much smaller than the order of the system.

Parameters

s

Complex number.

mu

Parameter values.

Returns

dtfs

Derivative of transfer function evaluated at the complex number s, NumPy array of shape (self.output_dim, self.input_dim).

eval_tf(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 either the number of inputs or the number of outputs is much smaller than the order of the system.

Parameters

s

Complex number.

mu

Parameter values.

Returns

tfs

Transfer function evaluated at the complex number s, NumPy array of shape (self.output_dim, self.input_dim).

classmethod from_files(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, 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.

estimator

An error estimator for the problem. This can be any object with an estimate(U, mu, model) method. If estimator is not None, an estimate(U, mu) method is added to the model which will call estimator.estimate(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

som

The SecondOrderModel with operators M, E, K, B, Cp, Cv, and D.

classmethod from_matrices(M, E, K, B, Cp, Cv=None, D=None, cont_time=True, state_id='STATE', solver_options=None, 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.

estimator

An error estimator for the problem. This can be any object with an estimate(U, mu, model) method. If estimator is not None, an estimate(U, mu) method is added to the model which will call estimator.estimate(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.

gramian(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.

h2_norm(mu=None)[source]

Compute the H2-norm.

Note

Assumes the system is asymptotically stable.

Parameters

mu

Parameter values.

Returns

norm

H_2-norm.

hankel_norm(mu=None)[source]

Compute the Hankel-norm.

Note

Assumes the system is asymptotically stable.

Parameters

mu

Parameter values.

Returns

norm

Hankel-norm.

hinf_norm(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).

poles(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.

psv(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.

pvsv(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.

to_lti()[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.

vpsv(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.

vsv(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.


class pymor.models.iosys.TransferFunction(*args, **kwargs)[source]

Bases: pymor.models.iosys.InputOutputModel

Class for systems represented by a transfer function.

This class describes input-output systems given by a transfer function \(H(s, mu)\).

Parameters

input_space

The input VectorSpace. Typically NumpyVectorSpace(m) where m is the number of inputs.

output_space

The output VectorSpace. Typically NumpyVectorSpace(p) where p is the number of outputs.

tf

The transfer function defined at least on the open right complex half-plane. tf(s, mu) is a NumPy array of shape (p, m).

dtf

The complex derivative of H with respect to s.

cont_time

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

name

Name of the system.

input_dim

The number of inputs.

output_dim

The number of outputs.

tf

The transfer function.

dtf

The complex derivative of the transfer function.

__str__()[source]

Return str(self).

eval_dtf(s, mu=None)[source]

Evaluate the derivative of the transfer function.

eval_tf(s, mu=None)[source]

Evaluate the transfer function.

h2_norm(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 argument epsabs). It can be changed by using the epsabs 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 is False).

info

Quadrature info (returned if return_norm_only is False and full_output is True). See scipy.integrate.quad documentation for more details.


pymor.models.iosys.sparse_min_size(value=1000)[source]

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

Defaults

value (see pymor.core.defaults)

mpi module


class pymor.models.mpi.MPIModel(obj_id, *args, **kwargs)[source]

Bases: object

Wrapper class mixin for MPI distributed Models.

Given a single-rank implementation of a Model, this wrapper class uses the event loop from pymor.tools.mpi to allow an MPI distributed usage of the Model. The underlying implementation needs to be MPI aware. In particular, the model’s solve method has to perform an MPI parallel solve of the model.

Note that this class is not intended to be instantiated directly. Instead, you should use mpi_wrap_model.

Methods

MPIModel

visualize


class pymor.models.mpi.MPIVisualizer(*args, **kwargs)[source]

Bases: pymor.core.base.ImmutableObject


class pymor.models.mpi._OperatorToWrap(operator, mpi_range, mpi_source)

Bases: tuple

Methods

_OperatorToWrap

__getnewargs__, __new__, __repr__

tuple

count, index

__getnewargs__()

Return self as a plain tuple. Used by copy and pickle.

static __new__(_cls, operator, mpi_range, mpi_source)

Create new instance of _OperatorToWrap(operator, mpi_range, mpi_source)

__repr__()

Return a nicely formatted representation string

property mpi_range

Alias for field number 1

property mpi_source

Alias for field number 2

property operator

Alias for field number 0


pymor.models.mpi.mpi_wrap_model(local_models, mpi_spaces=('STATE', ), use_with=True, with_apply2=False, pickle_local_spaces=True, space_type=<class 'pymor.vectorarrays.mpi.MPIVectorSpace'>, base_type=None)[source]

Wrap MPI distributed local Models to a global Model on rank 0.

Given MPI distributed local Models referred to by the ObjectId local_models, return a new Model which manages these distributed models from rank 0. This is done by first wrapping all Operators of the Model using mpi_wrap_operator.

Alternatively, local_models can be a callable (with no arguments) which is then called on each rank to instantiate the local Models.

When use_with is False, an MPIModel is instantiated with the wrapped operators. A call to solve will then use an MPI parallel call to the solve methods of the wrapped local Models to obtain the solution. This is usually what you want when the actual solve is performed by an implementation in the external solver.

When use_with is True, with_ is called on the local Model on rank 0, to obtain a new Model with the wrapped MPI Operators. This is mainly useful when the local models are generic Models as in pymor.models.basic and solve is implemented directly in pyMOR via operations on the contained Operators.

Parameters

local_models

ObjectId of the local Models on each rank or a callable generating the Models.

mpi_spaces

List of types or ids of VectorSpaces which are MPI distributed and need to be wrapped.

use_with

See above.

with_apply2

See MPIOperator.

pickle_local_spaces

See MPIOperator.

space_type

See MPIOperator.

neural_network module


class pymor.models.neural_network.FullyConnectedNN(layers_sizes, activation_function=<built-in method tanh of type object>)[source]

Bases: torch.nn.modules.module.Module, pymor.core.base.BasicObject

Class for neural networks with fully connected layers.

This class implements neural networks consisting of linear and fully connected layers. Furthermore, the same activation function is used between each layer, except for the last one where no activation function is applied.

Parameters

layers_sizes

List of sizes (i.e. number of neurons) for the layers of the neural network.

activation_function

Function to use as activation function between the single layers.

Methods

FullyConnectedNN

forward

Module

add_module, apply, bfloat16, buffers, children, cpu, cuda, double, eval, extra_repr, float, half, load_state_dict, modules, named_buffers, named_children, named_modules, named_parameters, parameters, register_backward_hook, register_buffer, register_forward_hook, register_forward_pre_hook, register_parameter, requires_grad_, share_memory, state_dict, to, train, type, zero_grad

BasicObject

disable_logging, enable_logging

Attributes

Module

dump_patches

BasicObject

logger, logging_disabled, name, uid

forward(x)[source]

Performs the forward pass through the neural network.

Applies the weights in the linear layers and passes the outcomes to the activation function.

Parameters

x

Input for the neural network.

Returns

The output of the neural network for the input x.


class pymor.models.neural_network.NeuralNetworkModel(*args, **kwargs)[source]

Bases: pymor.models.interface.Model

Class for models of stationary problems that use artificial neural networks.

This class implements a Model that uses a neural network for solving.

Parameters

neural_network

The neural network that approximates the mapping from parameter space to solution space. Should be an instance of FullyConnectedNN with input size that matches the (total) number of parameters and output size equal to the dimension of the reduced space.

output_functional

Operator mapping a given solution to the model output. In many applications, this will be a Functional, i.e. an Operator mapping to scalars. This is not required, however.

products

A dict of inner product Operators defined on the discrete space the problem is posed on. For each product with key 'x' a corresponding attribute x_product, as well as a norm method x_norm is added to the model.

estimator

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

visualizer

A visualizer for the problem. This can be any object with a visualize(U, m, ...) 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 model.