pymor.algorithms.lyapunov

Module Contents

Functions

mat_eqn_sparse_min_size

Returns minimal size for which a sparse solver will be used by default.

solve_lyap_lrcf

Compute an approximate low-rank solution of a Lyapunov equation.

_solve_lyap_lrcf_check_args

solve_lyap_dense

Compute the solution of a Lyapunov equation.

_solve_lyap_dense_check_args

_chol

Cholesky decomposition.

Attributes

_DEFAULT_LYAP_LRCF_SPARSE_SOLVER_BACKEND

_DEFAULT_LYAP_LRCF_DENSE_SOLVER_BACKEND

_DEFAULT_LYAP_DENSE_SOLVER_BACKEND

pymor.algorithms.lyapunov._DEFAULT_LYAP_LRCF_SPARSE_SOLVER_BACKEND[source]
pymor.algorithms.lyapunov._DEFAULT_LYAP_LRCF_DENSE_SOLVER_BACKEND[source]
pymor.algorithms.lyapunov._DEFAULT_LYAP_DENSE_SOLVER_BACKEND[source]
pymor.algorithms.lyapunov.mat_eqn_sparse_min_size(value=1000)[source]

Returns minimal size for which a sparse solver will be used by default.

pymor.algorithms.lyapunov.solve_lyap_lrcf(A, E, B, trans=False, options=None, default_sparse_solver_backend=_DEFAULT_LYAP_LRCF_SPARSE_SOLVER_BACKEND, default_dense_solver_backend=_DEFAULT_LYAP_LRCF_DENSE_SOLVER_BACKEND)[source]

Compute an approximate low-rank solution of a Lyapunov equation.

Returns a low-rank Cholesky factor \(Z\) such that \(Z Z^T\) approximates the solution \(X\) of a (generalized) continuous-time algebraic Lyapunov equation:

  • if trans is False and E is None:

    \[A X + X A^T + B B^T = 0,\]
  • if trans is False and E is an Operator:

    \[A X E^T + E X A^T + B B^T = 0,\]
  • if trans is True and E is None:

    \[A^T X + X A + B^T B = 0,\]
  • if trans is True and E is an Operator:

    \[A^T X E + E^T X A + B^T B = 0.\]

We assume A and E are real Operators, E is invertible, and all the eigenvalues of (A, E) all lie in the open left half-plane. Operator B needs to be given as a VectorArray from A.source, and for large-scale problems, we assume len(B) is small.

If the solver is not specified using the options argument, a solver backend is chosen based on availability in the following order:

Parameters

A

The non-parametric Operator A.

E

The non-parametric Operator E or None.

B

The operator B as a VectorArray from A.source.

trans

Whether the first Operator in the Lyapunov equation is transposed.

options

The solver options to use. See:

default_sparse_solver_backend

Default sparse solver backend to use (pymess, lradi).

default_dense_solver_backend

Default dense solver backend to use (pymess, slycot, scipy).

Returns

Z

Low-rank Cholesky factor of the Lyapunov equation solution, VectorArray from A.source.

pymor.algorithms.lyapunov._solve_lyap_lrcf_check_args(A, E, B, trans)[source]
pymor.algorithms.lyapunov.solve_lyap_dense(A, E, B, trans=False, options=None, default_solver_backend=_DEFAULT_LYAP_DENSE_SOLVER_BACKEND)[source]

Compute the solution of a Lyapunov equation.

Returns the solution \(X\) of a (generalized) continuous-time algebraic Lyapunov equation:

  • if trans is False and E is None:

    \[A X + X A^T + B B^T = 0,\]
  • if trans is False and E is an Operator:

    \[A X E^T + E X A^T + B B^T = 0,\]
  • if trans is True and E is None:

    \[A^T X + X A + B^T B = 0,\]
  • if trans is True and E is an Operator:

    \[A^T X E + E^T X A + B^T B = 0.\]

We assume A and E are real NumPy arrays, E is invertible, and that no two eigenvalues of (A, E) sum to zero (i.e., there exists a unique solution X).

If the solver is not specified using the options argument, a solver backend is chosen based on availability in the following order:

  1. pymess (see pymor.bindings.pymess.solve_lyap_dense)

  2. slycot (see pymor.bindings.slycot.solve_lyap_dense)

  3. scipy (see pymor.bindings.scipy.solve_lyap_dense)

Parameters

A

The matrix A as a 2D NumPy array.

E

The matrix E as a 2D NumPy array or None.

B

The matrix B as a 2D NumPy array.

trans

Whether the first operator in the Lyapunov equation is transposed.

options

The solver options to use. See:

default_solver_backend

Default solver backend to use (pymess, slycot, scipy).

Returns

X

Lyapunov equation solution as a NumPy array.

pymor.algorithms.lyapunov._solve_lyap_dense_check_args(A, E, B, trans)[source]
pymor.algorithms.lyapunov._chol(A)[source]

Cholesky decomposition.

This implementation uses SVD to compute the Cholesky factor (can be used for singular matrices).

Parameters

A

Symmetric positive semidefinite matrix as a NumPy array.

Returns

L

Cholesky factor of A (in the sense that L * L^T approximates A).