pymor.bindings.pymess¶
Module Contents¶
- class pymor.bindings.pymess.LyapunovEquation(opt, A, E, B)[source]¶
Bases:
pymess.EquationLyapunov equation class for pymess.
Represents a (generalized) continuous-time algebraic Lyapunov equation:
if opt.type is
pymess.MESS_OP_NONEand E isNone:\[A X + X A^T + B B^T = 0,\]if opt.type is
pymess.MESS_OP_NONEand E is notNone:\[A X E^T + E X A^T + B B^T = 0,\]if opt.type is
pymess.MESS_OP_TRANSPOSEand E isNone:\[A^T X + X A + B^T B = 0,\]if opt.type is
pymess.MESS_OP_TRANSPOSEand E is notNone:\[A^T X E + E^T X A + B^T B = 0.\]
Parameters
- opt
pymess Options structure.
- A
The non-parametric
OperatorA.- E
The non-parametric
OperatorE orNone.- B
The operator B as a
VectorArrayfromA.source.
- class pymor.bindings.pymess.RiccatiEquation(opt, A, E, B, C)[source]¶
Bases:
pymess.EquationRiccati equation class for pymess.
Represents a Riccati equation
if opt.type is
pymess.MESS_OP_NONEand E isNone:\[A X + X A^T - X C^T C X + B B^T = 0,\]if opt.type is
pymess.MESS_OP_NONEand E is notNone:\[A X E^T + E X A^T - E X C^T C X E^T + B B^T = 0,\]if opt.type is
pymess.MESS_OP_TRANSPOSEand E isNone:\[A^T X + X A - X B B^T X + C^T C = 0,\]if opt.type is
pymess.MESS_OP_TRANSPOSEand E is notNone:\[A^T X E + E^T X A - E X B B^T X E^T + C^T C = 0.\]
Parameters
- opt
pymess Options structure.
- 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.
- pymor.bindings.pymess.dense_nm_gmpcare_solver_options(linesearch=False, maxit=50, absres_tol=1e-11, relres_tol=1e-12, nrm=0)[source]¶
Return available Riccati solvers with default options for the pymess backend.
Also see
lradi_solver_options.Parameters
- linesearch
See
pymess.dense_nm_gmpcare.- maxit
See
pymess.dense_nm_gmpcare.- absres_tol
See
pymess.dense_nm_gmpcare.- relres_tol
See
pymess.dense_nm_gmpcare.- nrm
See
pymess.dense_nm_gmpcare.
Returns
A dict of available solvers with default solver options.
- pymor.bindings.pymess.lradi_solver_options(adi_maxit=500, adi_memory_usage=pymess.MESS_MEMORY_MID, adi_output=1, adi_rel_change_tol=1e-10, adi_res2_tol=1e-10, adi_res2c_tol=1e-11, adi_shifts_arp_m=32, adi_shifts_arp_p=48, adi_shifts_b0=None, adi_shifts_l0=16, adi_shifts_p=None, adi_shifts_paratype=pymess.MESS_LRCFADI_PARA_ADAPTIVE_Z)[source]¶
Return available adi solver options with default values for the pymess backend.
Parameters
- adi_maxit
See
pymess.OptionsAdi.- adi_memory_usage
See
pymess.OptionsAdi.- adi_output
See
pymess.OptionsAdi.- adi_rel_change_tol
See
pymess.OptionsAdi.- adi_res2_tol
See
pymess.OptionsAdi.- adi_res2c_tol
See
pymess.OptionsAdi.- adi_shifts_arp_m
See
pymess.OptionsAdiShifts.- adi_shifts_arp_p
See
pymess.OptionsAdiShifts.- adi_shifts_b0
See
pymess.OptionsAdiShifts.- adi_shifts_l0
See
pymess.OptionsAdiShifts.- adi_shifts_p
See
pymess.OptionsAdiShifts.- adi_shifts_paratype
See
pymess.OptionsAdiShifts.
Returns
A dict of available solvers with default solver options.
- pymor.bindings.pymess.lrnm_solver_options(newton_gstep=0, newton_k0=None, newton_linesearch=0, newton_maxit=30, newton_output=1, newton_res2_tol=1e-10, newton_singleshifts=0)[source]¶
Return available adi solver options with default values for the pymess backend.
Parameters
- newton_gstep
See
pymess.OptionsNewton.- newton_k0
See
pymess.OptionsNewton.- newton_linesearch
See
pymess.OptionsNewton.- newton_maxit
See
pymess.OptionsNewton.- newton_output
See
pymess.OptionsNewton.- newton_res2_tol
See
pymess.OptionsNewton.- newton_singleshifts
See
pymess.OptionsNewton.
Returns
A dict of available solvers with default solver options.
- pymor.bindings.pymess.lyap_dense_solver_options()[source]¶
Return available Lyapunov solvers with default options for the pymess backend.
Returns
A dict of available solvers with default solver options.
- pymor.bindings.pymess.lyap_lrcf_solver_options()[source]¶
Return available Lyapunov solvers with default options for the pymess backend.
Also see
lradi_solver_options.Returns
A dict of available solvers with default solver options.
- pymor.bindings.pymess.pos_ricc_lrcf_solver_options()[source]¶
Return available positive Riccati solvers with default options for the pymess backend.
Returns
A dict of available solvers with default solver options.
- pymor.bindings.pymess.ricc_lrcf_solver_options()[source]¶
Return available Riccati solvers with default options for the pymess backend.
Also see
dense_nm_gmpcare_solver_optionsandlrnm_solver_options.Returns
A dict of available solvers with default solver options.
- pymor.bindings.pymess.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
pymess.glyap.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.
- cont_time
Whether the continuous- or discrete-time Lyapunov equation is solved. Only the continuous-time case is implemented.
- options
The solver options to use (see
lyap_dense_solver_options).
Returns
- X
Lyapunov equation solution as a
NumPy array.
- pymor.bindings.pymess.solve_lyap_lrcf(A, E, B, trans=False, cont_time=True, options=None, default_solver=None)[source]¶
Compute an approximate low-rank solution of a Lyapunov equation.
See
for a general description.
This function uses
pymess.glyapandpymess.lradi. For both methods,to_numpyandfrom_numpyneed to be implemented forA.source. Additionally, sinceglyapis a dense solver, it expectsto_matrixto work for A and E.If the solver is not specified using the options or default_solver arguments,
glyapis used for small problems (smaller than defined withmat_eqn_sparse_min_size) andlradifor large problems.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.- cont_time
Whether the continuous- or discrete-time Lyapunov equation is solved. Only the continuous-time case is implemented.
- options
The solver options to use (see
lyap_lrcf_solver_options).- default_solver
Default solver to use (pymess_lradi, pymess_glyap). If
None, choose solver depending on the dimension of A.
Returns
- Z
Low-rank Cholesky factor of the Lyapunov equation solution,
VectorArrayfromA.source.
- pymor.bindings.pymess.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_lrcffor a general description.This function uses
pymess.dense_nm_gmpcare.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.- S
The operator S as a
VectorArrayfromA.sourceorNone.- trans
Whether the first
Operatorin the Riccati equation is transposed.- options
The solver options to use (see
pos_ricc_lrcf_solver_options).
Returns
- Z
Low-rank Cholesky factor of the Riccati equation solution,
VectorArrayfromA.source.
- pymor.bindings.pymess.solve_ricc_lrcf(A, E, B, C, R=None, S=None, trans=False, options=None, default_solver=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
pymess.dense_nm_gmpcareandpymess.lrnm. For both methods,to_numpyandfrom_numpyneed to be implemented forA.source. Additionally, sincedense_nm_gmpcareis a dense solver, it expectsto_matrixto work for A and E.If the solver is not specified using the options or default_solver arguments,
dense_nm_gmpcareis used for small problems (smaller than defined withmat_eqn_sparse_min_size) andlrnmfor large problems.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.- S
The operator S as a
VectorArrayfromA.sourceorNone.- trans
Whether the first
Operatorin the Riccati equation is transposed.- options
The solver options to use (see
ricc_lrcf_solver_options).- default_solver
Default solver to use (pymess_lrnm, pymess_dense_nm_gmpcare). If
None, chose solver depending on dimensionA.
Returns
- Z
Low-rank Cholesky factor of the Riccati equation solution,
VectorArrayfromA.source.