pymor.bindings.scipy

Module Contents

pymor.bindings.scipy.apply_inverse(op, V, initial_guess=None, options=None, least_squares=False, check_finite=True, default_solver='scipy_spsolve', default_least_squares_solver='scipy_least_squares_lsmr')[source]

Solve linear equation system.

Applies the inverse of op to the vectors in V using SciPy.

Parameters

op

The linear, non-parametric Operator to invert.

V

VectorArray of right-hand sides for the equation system.

initial_guess

VectorArray with the same length as V containing initial guesses for the solution. Some implementations of apply_inverse may ignore this parameter. If None a solver-dependent default is used.

options

The solver_options to use (see solver_options).

least_squares

If True, return least squares solution.

check_finite

Test if solution only contains finite values.

default_solver

Default solver to use (scipy_spsolve, scipy_bicgstab, scipy_bicgstab_spilu, scipy_lgmres, scipy_least_squares_lsmr, scipy_least_squares_lsqr).

default_least_squares_solver

Default solver to use for least squares problems (scipy_least_squares_lsmr, scipy_least_squares_lsqr).

Returns

VectorArray of the solution vectors.

pymor.bindings.scipy.lyap_dense_solver_options()[source]

Return available dense Lyapunov solvers with default options for the SciPy backend.

Returns

A dict of available solvers with default solver options.

pymor.bindings.scipy.lyap_lrcf_solver_options()[source]

Return available Lyapunov solvers with default options for the SciPy backend.

Returns

A dict of available solvers with default solver options.

pymor.bindings.scipy.matrix_astype_nocopy(matrix, dtype)[source]
pymor.bindings.scipy.ricc_dense_solver_options()[source]

Return available Riccati solvers with default options for the SciPy backend.

Returns

A dict of available solvers with default solver options.

pymor.bindings.scipy.ricc_lrcf_solver_options()[source]

Return available Riccati solvers with default options for the SciPy backend.

Returns

A dict of available solvers with default solver options.

pymor.bindings.scipy.solve_lyap_dense(A, E, B, trans=False, cont_time=True, options=None)[source]

Compute the solution of a Lyapunov equation.

See

for a general description.

This function uses scipy.linalg.solve_continuous_lyapunov or scipy.linalg.solve_discrete_lyapunov, which are dense solvers for Lyapunov equations with E=I.

Note

If E is not None, the problem will be reduced to a standard algebraic Lyapunov equation by inverting E.

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.

cont_time

Whether the continuous- or discrete-time Lyapunov equation is solved.

options

The solver options to use (see lyap_dense_solver_options).

Returns

X

Lyapunov equation solution as a NumPy array.

pymor.bindings.scipy.solve_lyap_lrcf(A, E, B, trans=False, cont_time=True, options=None)[source]

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

See

for a general description.

This function uses scipy.linalg.solve_continuous_lyapunov or scipy.linalg.solve_discrete_lyapunov, which are dense solvers for Lyapunov equations with E=I. Therefore, we assume A and E can be converted to NumPy arrays using to_matrix and that B.to_numpy is implemented.

Note

If E is not None, the problem will be reduced to a standard algebraic Lyapunov equation by inverting E.

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.

cont_time

Whether the continuous- or discrete-time Lyapunov equation is solved.

options

The solver options to use (see lyap_lrcf_solver_options).

Returns

Z

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

pymor.bindings.scipy.solve_pos_ricc_dense(A, E, B, C, R=None, S=None, trans=False, options=None)[source]

Compute the solution of a Riccati equation.

See pymor.algorithms.riccati.solve_pos_ricc_dense for a general description.

This function uses scipy.linalg.solve_continuous_are, which is a dense solver.

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.

C

The matrix C as a 2D NumPy array.

R

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

S

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

trans

Whether the first operator in the Riccati equation is transposed.

options

The solver options to use (see ricc_dense_solver_options).

Returns

X

Riccati equation solution as a NumPy array.

pymor.bindings.scipy.solve_pos_ricc_lrcf(A, E, B, C, R=None, S=None, trans=False, options=None)[source]

Compute an approximate low-rank solution of a positive Riccati equation.

See pymor.algorithms.riccati.solve_pos_ricc_lrcf for a general description.

This function uses scipy.linalg.solve_continuous_are, which is a dense solver. Therefore, we assume all Operators and VectorArrays can be converted to NumPy arrays using to_matrix and to_numpy.

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.

C

The operator C as a VectorArray from A.source.

R

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

S

The operator S as a VectorArray from A.source or None.

trans

Whether the first Operator in the positive Riccati equation is transposed.

options

The solver options to use (see ricc_lrcf_solver_options).

Returns

Z

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

pymor.bindings.scipy.solve_ricc_dense(A, E, B, C, R=None, S=None, trans=False, options=None)[source]

Compute the solution of a Riccati equation.

See pymor.algorithms.riccati.solve_ricc_dense for a general description.

This function uses scipy.linalg.solve_continuous_are, which is a dense solver.

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.

C

The matrix C as a 2D NumPy array.

R

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

S

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

trans

Whether the first operator in the Riccati equation is transposed.

options

The solver options to use (see ricc_dense_solver_options).

Returns

X

Riccati equation solution as a NumPy array.

pymor.bindings.scipy.solve_ricc_lrcf(A, E, B, C, R=None, S=None, trans=False, options=None)[source]

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

See pymor.algorithms.riccati.solve_ricc_lrcf for a general description.

This function uses scipy.linalg.solve_continuous_are, which is a dense solver. Therefore, we assume all Operators and VectorArrays can be converted to NumPy arrays using to_matrix and to_numpy.

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.

C

The operator C as a VectorArray from A.source.

R

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

S

The operator S as a VectorArray from A.source or None.

trans

Whether the first Operator in the Riccati equation is transposed.

options

The solver options to use (see ricc_lrcf_solver_options).

Returns

Z

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

pymor.bindings.scipy.solver_options(bicgstab_tol=1e-15, bicgstab_maxiter=None, spilu_drop_tol=0.0001, spilu_fill_factor=10, spilu_drop_rule=None, spilu_permc_spec='COLAMD', spsolve_permc_spec='COLAMD', spsolve_keep_factorization=True, lgmres_tol=1e-05, lgmres_maxiter=1000, lgmres_inner_m=39, lgmres_outer_k=3, least_squares_lsmr_damp=0.0, least_squares_lsmr_atol=1e-06, least_squares_lsmr_btol=1e-06, least_squares_lsmr_conlim=100000000.0, least_squares_lsmr_maxiter=None, least_squares_lsmr_show=False, least_squares_lsqr_damp=0.0, least_squares_lsqr_atol=1e-06, least_squares_lsqr_btol=1e-06, least_squares_lsqr_conlim=100000000.0, least_squares_lsqr_iter_lim=None, least_squares_lsqr_show=False)[source]

Returns available solvers with default solver_options for the SciPy backend.

Parameters

bicgstab_tol

See scipy.sparse.linalg.bicgstab.

bicgstab_maxiter

See scipy.sparse.linalg.bicgstab.

spilu_drop_tol

See scipy.sparse.linalg.spilu.

spilu_fill_factor

See scipy.sparse.linalg.spilu.

spilu_drop_rule

See scipy.sparse.linalg.spilu.

spilu_permc_spec

See scipy.sparse.linalg.spilu.

spsolve_permc_spec

See scipy.sparse.linalg.spsolve.

spsolve_keep_factorization

See scipy.sparse.linalg.spsolve.

lgmres_tol

See scipy.sparse.linalg.lgmres.

lgmres_maxiter

See scipy.sparse.linalg.lgmres.

lgmres_inner_m

See scipy.sparse.linalg.lgmres.

lgmres_outer_k

See scipy.sparse.linalg.lgmres.

least_squares_lsmr_damp

See scipy.sparse.linalg.lsmr.

least_squares_lsmr_atol

See scipy.sparse.linalg.lsmr.

least_squares_lsmr_btol

See scipy.sparse.linalg.lsmr.

least_squares_lsmr_conlim

See scipy.sparse.linalg.lsmr.

least_squares_lsmr_maxiter

See scipy.sparse.linalg.lsmr.

least_squares_lsmr_show

See scipy.sparse.linalg.lsmr.

least_squares_lsqr_damp

See scipy.sparse.linalg.lsqr.

least_squares_lsqr_atol

See scipy.sparse.linalg.lsqr.

least_squares_lsqr_btol

See scipy.sparse.linalg.lsqr.

least_squares_lsqr_conlim

See scipy.sparse.linalg.lsqr.

least_squares_lsqr_iter_lim

See scipy.sparse.linalg.lsqr.

least_squares_lsqr_show

See scipy.sparse.linalg.lsqr.

Returns

A dict of available solvers with default solver_options.

pymor.bindings.scipy.svd_lapack_driver(driver='gesvd_unless_win_mkl')[source]
pymor.bindings.scipy.SCIPY_1_14_OR_NEWER[source]
pymor.bindings.scipy.sparray[source]