pymor.bindings.scipy¶
Module Contents¶
Functions¶
Returns available solvers with default |
|
Solve linear equation system. |
|
Return available Lyapunov solvers with default options for the SciPy backend. |
|
Compute an approximate low-rank solution of a Lyapunov equation. |
|
Return available dense Lyapunov solvers with default options for the SciPy backend. |
|
Compute the solution of a Lyapunov equation. |
|
Return available Riccati solvers with default options for the SciPy backend. |
|
Compute an approximate low-rank solution of a Riccati equation. |
|
Return available Riccati solvers with default options for the SciPy backend. |
|
Compute the solution of a Riccati equation. |
|
Return available positive Riccati solvers with default options for the SciPy backend. |
|
Compute an approximate low-rank solution of a positive Riccati equation. |
- 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_optionsfor the SciPy backend.Parameters
- bicgstab_tol
- bicgstab_maxiter
- spilu_drop_tol
- spilu_fill_factor
- spilu_drop_rule
- spilu_permc_spec
- spsolve_permc_spec
- spsolve_keep_factorization
- lgmres_tol
- lgmres_maxiter
- lgmres_inner_m
- lgmres_outer_k
- least_squares_lsmr_damp
- least_squares_lsmr_atol
- least_squares_lsmr_btol
- least_squares_lsmr_conlim
- least_squares_lsmr_maxiter
- least_squares_lsmr_show
- least_squares_lsqr_damp
- least_squares_lsqr_atol
- least_squares_lsqr_btol
- least_squares_lsqr_conlim
- least_squares_lsqr_iter_lim
- least_squares_lsqr_show
Returns
A dict of available solvers with default
solver_options.
- 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
opto the vectors inVusing SciPy.Parameters
- op
The linear, non-parametric
Operatorto invert.- V
VectorArrayof right-hand sides for the equation system.- initial_guess
VectorArraywith the same length asVcontaining initial guesses for the solution. Some implementations ofapply_inversemay ignore this parameter. IfNonea solver-dependent default is used.- options
The
solver_optionsto use (seesolver_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
VectorArrayof the solution vectors.
- 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.solve_lyap_lrcf(A, E, B, trans=False, options=None)[source]¶
Compute an approximate low-rank solution of a Lyapunov equation.
See
pymor.algorithms.lyapunov.solve_lyap_lrcffor a general description.This function uses
scipy.linalg.solve_continuous_lyapunov, which is a dense solver for Lyapunov equations with E=I. Therefore, we assume A and E can be converted toNumPy arraysusingto_matrixand thatB.to_numpyis implemented.Note
If E is not
None, the problem will be reduced to a standard continuous-time algebraic Lyapunov equation by inverting E.Parameters
- A
The non-parametric
OperatorA.- E
The non-parametric
OperatorE orNone.- B
The operator B as a
VectorArrayfromA.source.- trans
Whether the first
Operatorin the Lyapunov equation is transposed.- options
The solver options to use (see
lyap_lrcf_solver_options).
Returns
- Z
Low-rank Cholesky factor of the Lyapunov equation solution,
VectorArrayfromA.source.
- 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.solve_lyap_dense(A, E, B, trans=False, options=None)[source]¶
Compute the solution of a Lyapunov equation.
See
pymor.algorithms.lyapunov.solve_lyap_densefor a general description.This function uses
scipy.linalg.solve_continuous_lyapunov, which is a dense solver for Lyapunov equations with E=I.Note
If E is not
None, the problem will be reduced to a standard continuous-time algebraic Lyapunov equation by inverting E.Parameters
- A
The matrix A as a 2D
NumPy array.- E
The matrix E as a 2D
NumPy arrayorNone.- 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
lyap_dense_solver_options).
Returns
- X
Lyapunov equation solution as a
NumPy array.
- 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_ricc_lrcf(A, E, B, C, R=None, trans=False, options=None)[source]¶
Compute an approximate low-rank solution of a Riccati equation.
See
pymor.algorithms.riccati.solve_ricc_lrcffor a general description.This function uses
scipy.linalg.solve_continuous_are, which is a dense solver. Therefore, we assume allOperatorsandVectorArrayscan be converted toNumPy arraysusingto_matrixandto_numpy.Parameters
- A
The non-parametric
OperatorA.- E
The non-parametric
OperatorE orNone.- B
The operator B as a
VectorArrayfromA.source.- C
The operator C as a
VectorArrayfromA.source.- R
The matrix R as a 2D
NumPy arrayorNone.- trans
Whether the first
Operatorin 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,
VectorArrayfromA.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.solve_ricc_dense(A, E, B, C, R=None, trans=False, options=None)[source]¶
Compute the solution of a Riccati equation.
See
pymor.algorithms.riccati.solve_ricc_densefor 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 arrayorNone.- 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 arrayorNone.- 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.pos_ricc_lrcf_solver_options()[source]¶
Return available positive Riccati solvers with default options for the SciPy backend.
Returns
A dict of available solvers with default solver options.
- pymor.bindings.scipy.solve_pos_ricc_lrcf(A, E, B, C, R=None, trans=False, options=None)[source]¶
Compute an approximate low-rank solution of a positive Riccati equation.
See
pymor.algorithms.riccati.solve_pos_ricc_lrcffor a general description.This function uses
scipy.linalg.solve_continuous_are, which is a dense solver. Therefore, we assume allOperatorsandVectorArrayscan be converted toNumPy arraysusingto_matrixandto_numpy.Parameters
- A
The non-parametric
OperatorA.- E
The non-parametric
OperatorE orNone.- B
The operator B as a
VectorArrayfromA.source.- C
The operator C as a
VectorArrayfromA.source.- R
The matrix R as a 2D
NumPy arrayorNone.- trans
Whether the first
Operatorin the positive Riccati equation is transposed.- options
The solver options to use (see
pos_ricc_lrcf_solver_options).
Returns
- Z
Low-rank Cholesky factor of the positive Riccati equation solution,
VectorArrayfromA.source.