pymor.bindings package¶
Submodules¶
fenics module¶
-
class
pymor.bindings.fenics.
ComplexifiedFenicsVector
(real_part, imag_part)[source]¶ Bases:
pymor.vectorarrays.list.ComplexifiedVector
Methods
amax
axpy
,conj
,copy
,dofs
,dot
,l1_norm
,l2_norm
,l2_norm2
,scal
,sup_norm
,to_numpy
Attributes
imag
,real
-
class
pymor.bindings.fenics.
FenicsMatrixOperator
(*args, **kwargs)[source]¶ Bases:
pymor.operators.list.LinearComplexifiedListVectorArrayOperatorBase
Wraps a FEniCS matrix as an
Operator
.Methods
apply2
,as_range_array
,as_source_array
,as_vector
,assemble
,d_mu
,jacobian
,pairwise_apply2
,restricted
,__add__
,__matmul__
,__mul__
,__radd__
-
class
pymor.bindings.fenics.
FenicsOperator
(*args, **kwargs)[source]¶ Bases:
pymor.operators.interface.Operator
Wraps a FEniCS form as an
Operator
.Methods
apply2
,apply_adjoint
,apply_inverse
,apply_inverse_adjoint
,as_range_array
,as_source_array
,as_vector
,assemble
,d_mu
,pairwise_apply2
,__add__
,__matmul__
,__mul__
,__radd__
Attributes
parameters
,parameters_inherited
,parameters_internal
,parameters_own
,parametric
-
apply
(U, mu=None)[source]¶ Apply the operator to a
VectorArray
.Parameters
- U
VectorArray
of vectors to which the operator is applied.- mu
The
parameter values
for which to evaluate the operator.
Returns
VectorArray
of the operator evaluations.
-
jacobian
(U, mu=None)[source]¶ Return the operator’s Jacobian as a new
Operator
.Parameters
- U
Length 1
VectorArray
containing the vector for which to compute the Jacobian.- mu
The
parameter values
for which to compute the Jacobian.
Returns
Linear
Operator
representing the Jacobian.
-
restricted
(dofs)[source]¶ Restrict the operator range to a given set of degrees of freedom.
This method returns a restricted version
restricted_op
of the operator along with an arraysource_dofs
such that for anyVectorArray
U
inself.source
the following is true:self.apply(U, mu).dofs(dofs) == restricted_op.apply(NumpyVectorArray(U.dofs(source_dofs)), mu))
Such an operator is mainly useful for
empirical interpolation
where the evaluation of the original operator only needs to be known for few selected degrees of freedom. If the operator has a small stencil, only fewsource_dofs
will be needed to evaluate the restricted operator which can make its evaluation very fast compared to evaluating the original operator.Parameters
- dofs
One-dimensional
NumPy array
of degrees of freedom in the operatorrange
to which to restrict.
Returns
- restricted_op
The restricted operator as defined above. The operator will have
NumpyVectorSpace
(len(source_dofs))
assource
andNumpyVectorSpace
(len(dofs))
asrange
.- source_dofs
One-dimensional
NumPy array
of source degrees of freedom as defined above.
-
-
class
pymor.bindings.fenics.
FenicsVector
(impl)[source]¶ Bases:
pymor.vectorarrays.list.CopyOnWriteVector
Wraps a FEniCS vector to make it usable with ListVectorArray.
Methods
amax
,dofs
,dot
,from_instance
,l1_norm
,l2_norm
,l2_norm2
,sup_norm
,to_numpy
axpy
,copy
,scal
conj
Attributes
imag
,real
-
class
pymor.bindings.fenics.
FenicsVectorSpace
(*args, **kwargs)[source]¶ Bases:
pymor.vectorarrays.list.ComplexifiedListVectorSpace
Methods
real_full_vector
,real_make_vector
,real_random_vector
,real_vector_from_numpy
,real_zero_vector
full_vector
,make_vector
,random_vector
,vector_from_numpy
,zero_vector
from_numpy
,full
,make_array
,ones
,ones_vector
,random
,space_from_dim
,space_from_vector_obj
,zeros
empty
,from_data
Attributes
-
complexified_vector_type
¶ alias of
ComplexifiedFenicsVector
-
-
class
pymor.bindings.fenics.
FenicsVisualizer
(space, mesh_refinements=0)[source]¶ Bases:
pymor.core.base.BasicObject
Visualize a FEniCS grid function.
Parameters
- space
The
FenicsVectorSpace
for which we want to visualize DOF vectors.- mesh_refinements
Number of uniform mesh refinements to perform for vtk visualization (of functions from higher-order FE spaces).
Attributes
-
visualize
(U, m, title='', legend=None, filename=None, block=True, separate_colorbars=True)[source]¶ Visualize the provided data.
Parameters
- U
VectorArray
of the data to visualize (length must be 1). Alternatively, a tuple ofVectorArrays
which will be visualized in separate windows. Iffilename
is specified, only oneVectorArray
may be provided which, however, is allowed to contain multipled vectors that will be interpreted as a time series.- m
Filled in by
pymor.models.interface.Model.visualize
(ignored).- title
Title of the plot.
- legend
Description of the data that is plotted. If
U
is a tuple ofVectorArrays
,legend
has to be a tuple of the same length.- filename
If specified, write the data to that file.
filename
needs to have an extension supported by FEniCS (e.g..pvd
).- separate_colorbars
If
True
, use separate colorbars for each subplot.- block
If
True
, block execution until the plot window is closed.
-
class
pymor.bindings.fenics.
RestrictedFenicsOperator
(*args, **kwargs)[source]¶ Bases:
pymor.operators.interface.Operator
Methods
apply2
,apply_adjoint
,apply_inverse
,apply_inverse_adjoint
,as_range_array
,as_source_array
,as_vector
,assemble
,d_mu
,pairwise_apply2
,restricted
,__add__
,__matmul__
,__mul__
,__radd__
Attributes
parameters
,parameters_inherited
,parameters_internal
,parameters_own
,parametric
-
apply
(U, mu=None)[source]¶ Apply the operator to a
VectorArray
.Parameters
- U
VectorArray
of vectors to which the operator is applied.- mu
The
parameter values
for which to evaluate the operator.
Returns
VectorArray
of the operator evaluations.
-
jacobian
(U, mu=None)[source]¶ Return the operator’s Jacobian as a new
Operator
.Parameters
- U
Length 1
VectorArray
containing the vector for which to compute the Jacobian.- mu
The
parameter values
for which to compute the Jacobian.
Returns
Linear
Operator
representing the Jacobian.
-
ngsolve module¶
-
class
pymor.bindings.ngsolve.
ComplexifiedNGSolveVector
(real_part, imag_part)[source]¶ Bases:
pymor.bindings.ngsolve.NGSolveVectorCommon
,pymor.vectorarrays.list.ComplexifiedVector
Methods
amax
,dofs
,l1_norm
axpy
,conj
,copy
,dot
,l2_norm
,l2_norm2
,scal
,sup_norm
,to_numpy
Attributes
imag
,real
-
class
pymor.bindings.ngsolve.
NGSolveMatrixOperator
(*args, **kwargs)[source]¶ Bases:
pymor.operators.list.LinearComplexifiedListVectorArrayOperatorBase
Wraps a NGSolve matrix as an
Operator
.Methods
apply2
,as_range_array
,as_source_array
,assemble
,d_mu
,jacobian
,pairwise_apply2
,restricted
,__add__
,__matmul__
,__mul__
,__radd__
Attributes
parameters
,parameters_inherited
,parameters_internal
,parameters_own
,parametric
-
as_vector
(copy=True)[source]¶ Return a vector representation of a linear functional or vector operator.
Depending on the operator’s
source
andrange
, this method is equivalent to callingas_range_array
oras_source_array
respectively. The resultingVectorArray
is required to have length 1.Parameters
- mu
The
parameter values
for which to return the vector representation.
Returns
- V
VectorArray
of length 1 containing the vector representation.
-
-
class
pymor.bindings.ngsolve.
NGSolveVector
(impl)[source]¶ Bases:
pymor.bindings.ngsolve.NGSolveVectorCommon
,pymor.vectorarrays.list.CopyOnWriteVector
Wraps a NGSolve BaseVector to make it usable with ListVectorArray.
Methods
dot
,from_instance
,l2_norm
,l2_norm2
,to_numpy
amax
,dofs
,l1_norm
axpy
,copy
,scal
conj
,sup_norm
Attributes
imag
,real
-
class
pymor.bindings.ngsolve.
NGSolveVectorSpace
(*args, **kwargs)[source]¶ Bases:
pymor.vectorarrays.list.ComplexifiedListVectorSpace
Methods
real_make_vector
,real_vector_from_numpy
,real_zero_vector
,space_from_vector_obj
full_vector
,make_vector
,random_vector
,real_full_vector
,real_random_vector
,vector_from_numpy
,zero_vector
from_numpy
,full
,make_array
,ones
,ones_vector
,random
,space_from_dim
,zeros
empty
,from_data
Attributes
complexified_vector_type
,dim
,value_dim
-
complexified_vector_type
¶ alias of
ComplexifiedNGSolveVector
-
-
class
pymor.bindings.ngsolve.
NGSolveVisualizer
(*args, **kwargs)[source]¶ Bases:
pymor.core.base.ImmutableObject
Visualize an NGSolve grid function.
Methods
Attributes
pyamg module¶
-
pymor.bindings.pyamg.
apply_inverse
(op, V, initial_guess=None, options=None, least_squares=False, check_finite=True, default_solver='pyamg_solve')[source]¶ Solve linear equation system.
Applies the inverse of
op
to the vectors inV
using PyAMG.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 asV
containing initial guesses for the solution. Some implementations ofapply_inverse
may ignore this parameter. IfNone
a solver-dependent default is used.- options
The
solver_options
to use (seesolver_options
).- least_squares
Must be
False
.- check_finite
Test if solution only contains finite values.
- default_solver
Default solver to use (pyamg_solve, pyamg_rs, pyamg_sa).
Returns
VectorArray
of the solution vectors.Defaults
check_finite, default_solver (see
pymor.core.defaults
)
-
pymor.bindings.pyamg.
solver_options
(tol=1e-05, maxiter=400, verb=False, rs_strength='classical', {'theta': 0.25}, rs_CF='RS', rs_presmoother='gauss_seidel', {'sweep': 'symmetric'}, rs_postsmoother='gauss_seidel', {'sweep': 'symmetric'}, rs_max_levels=10, rs_max_coarse=500, rs_coarse_solver='pinv2', rs_cycle='V', rs_accel=None, rs_tol=1e-05, rs_maxiter=100, sa_symmetry='hermitian', sa_strength='symmetric', sa_aggregate='standard', sa_smooth='jacobi', {'omega': 1.3333333333333333}, sa_presmoother='block_gauss_seidel', {'sweep': 'symmetric'}, sa_postsmoother='block_gauss_seidel', {'sweep': 'symmetric'}, sa_improve_candidates='block_gauss_seidel', {'sweep': 'symmetric', 'iterations': 4}, None, sa_max_levels=10, sa_max_coarse=500, sa_diagonal_dominance=False, sa_coarse_solver='pinv2', sa_cycle='V', sa_accel=None, sa_tol=1e-05, sa_maxiter=100)[source]¶ Returns available solvers with default
solver_options
for the PyAMG backend.Parameters
- tol
Tolerance for PyAMG blackbox solver.
- maxiter
Maximum iterations for PyAMG blackbox solver.
- verb
Verbosity flag for PyAMG blackbox solver.
- rs_strength
Parameter for PyAMG Ruge-Stuben solver.
- rs_CF
Parameter for PyAMG Ruge-Stuben solver.
- rs_presmoother
Parameter for PyAMG Ruge-Stuben solver.
- rs_postsmoother
Parameter for PyAMG Ruge-Stuben solver.
- rs_max_levels
Parameter for PyAMG Ruge-Stuben solver.
- rs_max_coarse
Parameter for PyAMG Ruge-Stuben solver.
- rs_coarse_solver
Parameter for PyAMG Ruge-Stuben solver.
- rs_cycle
Parameter for PyAMG Ruge-Stuben solver.
- rs_accel
Parameter for PyAMG Ruge-Stuben solver.
- rs_tol
Parameter for PyAMG Ruge-Stuben solver.
- rs_maxiter
Parameter for PyAMG Ruge-Stuben solver.
- sa_symmetry
Parameter for PyAMG Smoothed-Aggregation solver.
- sa_strength
Parameter for PyAMG Smoothed-Aggregation solver.
- sa_aggregate
Parameter for PyAMG Smoothed-Aggregation solver.
- sa_smooth
Parameter for PyAMG Smoothed-Aggregation solver.
- sa_presmoother
Parameter for PyAMG Smoothed-Aggregation solver.
- sa_postsmoother
Parameter for PyAMG Smoothed-Aggregation solver.
- sa_improve_candidates
Parameter for PyAMG Smoothed-Aggregation solver.
- sa_max_levels
Parameter for PyAMG Smoothed-Aggregation solver.
- sa_max_coarse
Parameter for PyAMG Smoothed-Aggregation solver.
- sa_diagonal_dominance
Parameter for PyAMG Smoothed-Aggregation solver.
- sa_coarse_solver
Parameter for PyAMG Smoothed-Aggregation solver.
- sa_cycle
Parameter for PyAMG Smoothed-Aggregation solver.
- sa_accel
Parameter for PyAMG Smoothed-Aggregation solver.
- sa_tol
Parameter for PyAMG Smoothed-Aggregation solver.
- sa_maxiter
Parameter for PyAMG Smoothed-Aggregation solver.
Returns
A dict of available solvers with default
solver_options
.Defaults
tol, maxiter, verb, rs_strength, rs_CF, rs_postsmoother, rs_max_levels, rs_max_coarse, rs_coarse_solver, rs_cycle, rs_accel, rs_tol, rs_maxiter, sa_symmetry, sa_strength, sa_aggregate, sa_smooth, sa_presmoother, sa_postsmoother, sa_improve_candidates, sa_max_levels, sa_max_coarse, sa_diagonal_dominance, sa_coarse_solver, sa_cycle, sa_accel, sa_tol, sa_maxiter (see
pymor.core.defaults
)
pymess module¶
-
class
pymor.bindings.pymess.
LyapunovEquation
(opt, A, E, B)[source]¶ Bases:
pymess.equation.Equation
Lyapunov equation class for pymess
Represents a (generalized) continuous-time algebraic Lyapunov equation:
if opt.type is
pymess.MESS_OP_NONE
and E isNone
:\[A X + X A^T + B B^T = 0,\]if opt.type is
pymess.MESS_OP_NONE
and E is notNone
:\[A X E^T + E X A^T + B B^T = 0,\]if opt.type is
pymess.MESS_OP_TRANSPOSE
and E isNone
:\[A^T X + X A + B^T B = 0,\]if opt.type is
pymess.MESS_OP_TRANSPOSE
and E is notNone
:\[A^T X E + E^T X A + B^T B = 0.\]
Parameters
- opt
pymess Options structure.
- A
The non-parametric
Operator
A.- E
The non-parametric
Operator
E orNone
.- B
The operator B as a
VectorArray
fromA.source
.
Methods
ainv_apply
,apeinv_apply
,apex_apply
,ax_apply
,einv_apply
,ex_apply
,parameter
Equation
ainv_clear
,ainv_generate
,apeinv_clear
,apeinv_generate
,apex_clear
,apex_generate
,ax_clear
,ax_generate
,einv_clear
,einv_generate
,ex_clear
,ex_generate
,__repr__
,__str__
Attributes
Equation
b
,c
,dim
,k
,name
,options
,rhs
-
ainv_apply
(op, y)[source]¶ Applies the inverse of \(A\) to a right hand side. It has to return the solution of \(Ax =y\).
-
apeinv_apply
(op, p, idx_p, y)[source]¶ Applies the inverse of \(A+pE\) to a right hand side. It has to return the solution of \((A+pE)x = y\).
-
apex_apply
(op, p, idx_p, y)[source]¶ Applies function the operator of \(A+pE\). It has to return \(x=(A+pE)y\).
-
ax_apply
(op, y)[source]¶ Applies \(A\) to a right hand side. It has to return the result of \(x = Ay\).
-
einv_apply
(op, y)[source]¶ Applies the inverse of \(E\) to a right hand side. It has to return the solution of \(Ex = y\).
-
class
pymor.bindings.pymess.
RiccatiEquation
(opt, A, E, B, C)[source]¶ Bases:
pymess.equation.Equation
Riccati equation class for pymess
Represents a Riccati equation
if opt.type is
pymess.MESS_OP_NONE
and E isNone
:\[A X + X A^T - X C^T C X + B B^T = 0,\]if opt.type is
pymess.MESS_OP_NONE
and 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_TRANSPOSE
and E isNone
:\[A^T X + X A - X B B^T X + C^T C = 0,\]if opt.type is
pymess.MESS_OP_TRANSPOSE
and 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
Operator
A.- E
The non-parametric
Operator
E orNone
.- B
The operator B as a
VectorArray
fromA.source
.- C
The operator C as a
VectorArray
fromA.source
.
Methods
ainv_apply
,apeinv_apply
,apex_apply
,ax_apply
,einv_apply
,ex_apply
,parameter
Equation
ainv_clear
,ainv_generate
,apeinv_clear
,apeinv_generate
,apex_clear
,apex_generate
,ax_clear
,ax_generate
,einv_clear
,einv_generate
,ex_clear
,ex_generate
,__repr__
,__str__
Attributes
Equation
b
,c
,dim
,k
,name
,options
,rhs
-
ainv_apply
(op, y)[source]¶ Applies the inverse of \(A\) to a right hand side. It has to return the solution of \(Ax =y\).
-
apeinv_apply
(op, p, idx_p, y)[source]¶ Applies the inverse of \(A+pE\) to a right hand side. It has to return the solution of \((A+pE)x = y\).
-
apex_apply
(op, p, idx_p, y)[source]¶ Applies function the operator of \(A+pE\). It has to return \(x=(A+pE)y\).
-
ax_apply
(op, y)[source]¶ Applies \(A\) to a right hand side. It has to return the result of \(x = Ay\).
-
einv_apply
(op, y)[source]¶ Applies the inverse of \(E\) to a right hand side. It has to return the solution of \(Ex = y\).
-
pymor.bindings.pymess.
_call_pymess_dense_nm_gmpare
(A, E, B, C, R, S, trans=False, options=None, plus=False, method_name='')[source]¶ Return the solution from pymess.dense_nm_gmpare solver.
-
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.
Defaults
linesearch, maxit, absres_tol, relres_tol, nrm (see
pymor.core.defaults
)
-
pymor.bindings.pymess.
lradi_solver_options
(adi_maxit=500, adi_memory_usage=1, 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=3)[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.
Defaults
adi_maxit, adi_memory_usage, adi_output, adi_rel_change_tol, adi_res2_tol, adi_res2c_tol, adi_shifts_arp_m, adi_shifts_arp_p, adi_shifts_b0, adi_shifts_l0, adi_shifts_p, adi_shifts_paratype (see
pymor.core.defaults
)
-
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.
Defaults
newton_gstep, newton_k0, newton_linesearch, newton_maxit, newton_output, newton_res2_tol, newton_singleshifts (see
pymor.core.defaults
)
-
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_options
andlrnm_solver_options
.Returns
A dict of available solvers with default solver options.
-
pymor.bindings.pymess.
solve_lyap_dense
(A, E, B, trans=False, options=None)[source]¶ Compute the solution of a Lyapunov equation.
See
pymor.algorithms.lyapunov.solve_lyap_dense
for a general description.This function uses
pymess.glyap
.Parameters
- A
The operator A as a 2D
NumPy array
.- E
The operator E as a 2D
NumPy array
orNone
.- B
The operator 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.pymess.
solve_lyap_lrcf
(A, E, B, trans=False, options=None, default_solver=None)[source]¶ Compute an approximate low-rank solution of a Lyapunov equation.
See
pymor.algorithms.lyapunov.solve_lyap_lrcf
for a general description.This function uses
pymess.glyap
andpymess.lradi
. For both methods,to_numpy
andfrom_numpy
need to be implemented forA.source
. Additionally, sinceglyap
is a dense solver, it expectsto_matrix
to work for A and E.If the solver is not specified using the options or default_solver arguments,
glyap
is used for small problems (smaller than defined withmat_eqn_sparse_min_size
) andlradi
for large problems.Parameters
- A
The non-parametric
Operator
A.- E
The non-parametric
Operator
E orNone
.- B
The operator B as a
VectorArray
fromA.source
.- trans
Whether the first
Operator
in the Lyapunov equation is transposed.- 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,
VectorArray
fromA.source
.
Defaults
default_solver (see
pymor.core.defaults
)
-
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_lrcf
for a general description.This function uses
pymess.dense_nm_gmpcare
.Parameters
- A
The non-parametric
Operator
A.- E
The non-parametric
Operator
E orNone
.- B
The operator B as a
VectorArray
fromA.source
.- C
The operator C as a
VectorArray
fromA.source
.- R
The operator R as a 2D
NumPy array
orNone
.- S
The operator S as a
VectorArray
fromA.source
orNone
.- trans
Whether the first
Operator
in 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,
VectorArray
fromA.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_lrcf
for a general description.This function uses
pymess.dense_nm_gmpcare
andpymess.lrnm
. For both methods,to_numpy
andfrom_numpy
need to be implemented forA.source
. Additionally, sincedense_nm_gmpcare
is a dense solver, it expectsto_matrix
to work for A and E.If the solver is not specified using the options or default_solver arguments,
dense_nm_gmpcare
is used for small problems (smaller than defined withmat_eqn_sparse_min_size
) andlrnm
for large problems.Parameters
- A
The non-parametric
Operator
A.- E
The non-parametric
Operator
E orNone
.- B
The operator B as a
VectorArray
fromA.source
.- C
The operator C as a
VectorArray
fromA.source
.- R
The operator R as a 2D
NumPy array
orNone
.- S
The operator S as a
VectorArray
fromA.source
orNone
.- trans
Whether the first
Operator
in 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,
VectorArray
fromA.source
.
Defaults
default_solver (see
pymor.core.defaults
)
scipy module¶
-
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 inV
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 asV
containing initial guesses for the solution. Some implementations ofapply_inverse
may ignore this parameter. IfNone
a solver-dependent default is used.- options
The
solver_options
to 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
VectorArray
of the solution vectors.Defaults
check_finite, default_solver, default_least_squares_solver (see
pymor.core.defaults
)
-
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.
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.
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, options=None)[source]¶ Compute the solution of a Lyapunov equation.
See
pymor.algorithms.lyapunov.solve_lyap_dense
for 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 operator A as a 2D
NumPy array
.- E
The operator E as a 2D
NumPy array
orNone
.- B
The operator 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.
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_lrcf
for 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 arrays
usingto_matrix
and thatB.to_numpy
is 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
Operator
A.- E
The non-parametric
Operator
E orNone
.- B
The operator B as a
VectorArray
fromA.source
.- trans
Whether the first
Operator
in 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,
VectorArray
fromA.source
.
-
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 allOperators
andVectorArrays
can be converted toNumPy arrays
usingto_matrix
andto_numpy
.Parameters
- A
The non-parametric
Operator
A.- E
The non-parametric
Operator
E orNone
.- B
The operator B as a
VectorArray
fromA.source
.- C
The operator C as a
VectorArray
fromA.source
.- R
The operator R as a 2D
NumPy array
orNone
.- S
The operator S as a
VectorArray
fromA.source
orNone
.- trans
Whether the first
Operator
in 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,
VectorArray
fromA.source
.
-
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 allOperators
andVectorArrays
can be converted toNumPy arrays
usingto_matrix
andto_numpy
.Parameters
- A
The non-parametric
Operator
A.- E
The non-parametric
Operator
E orNone
.- B
The operator B as a
VectorArray
fromA.source
.- C
The operator C as a
VectorArray
fromA.source
.- R
The operator R as a 2D
NumPy array
orNone
.- S
The operator S as a
VectorArray
fromA.source
orNone
.- 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
fromA.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
- 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
.Defaults
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_atol, least_squares_lsqr_btol, least_squares_lsqr_conlim, least_squares_lsqr_iter_lim, least_squares_lsqr_show (see
pymor.core.defaults
)
slycot module¶
-
pymor.bindings.slycot.
lyap_dense_solver_options
()[source]¶ Return available Lyapunov solvers with default options for the slycot backend.
Returns
A dict of available solvers with default solver options.
-
pymor.bindings.slycot.
lyap_lrcf_solver_options
()[source]¶ Return available Lyapunov solvers with default options for the slycot backend.
Returns
A dict of available solvers with default solver options.
-
pymor.bindings.slycot.
pos_ricc_lrcf_solver_options
()[source]¶ Return available positive Riccati solvers with default options for the slycot backend.
Returns
A dict of available solvers with default solver options.
-
pymor.bindings.slycot.
ricc_lrcf_solver_options
()[source]¶ Return available Riccati solvers with default options for the slycot backend.
Returns
A dict of available solvers with default solver options.
-
pymor.bindings.slycot.
solve_lyap_dense
(A, E, B, trans=False, options=None)[source]¶ Compute the solution of a Lyapunov equation.
See
pymor.algorithms.lyapunov.solve_lyap_dense
for a general description.This function uses
slycot.sb03md
(ifE is None
) andslycot.sg03ad
(ifE is not None
), which are based on the Bartels-Stewart algorithm.Parameters
- A
The operator A as a 2D
NumPy array
.- E
The operator E as a 2D
NumPy array
orNone
.- B
The operator 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.slycot.
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_lrcf
for a general description.This function uses
slycot.sb03md
(ifE is None
) andslycot.sg03ad
(ifE is not None
), which are dense solvers based on the Bartels-Stewart algorithm. Therefore, we assume A and E can be converted toNumPy arrays
usingto_matrix
and thatB.to_numpy
is implemented.Parameters
- A
The non-parametric
Operator
A.- E
The non-parametric
Operator
E orNone
.- B
The operator B as a
VectorArray
fromA.source
.- trans
Whether the first
Operator
in 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,
VectorArray
fromA.source
.
-
pymor.bindings.slycot.
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
slycot.sb02md
(if E and S areNone
),slycot.sb02od
(if E isNone
and S is notNone
) andslycot.sg03ad
(if E is notNone
), which are dense solvers. Therefore, we assume allOperators
andVectorArrays
can be converted toNumPy arrays
usingto_matrix
andto_numpy
.Parameters
- A
The non-parametric
Operator
A.- E
The non-parametric
Operator
E orNone
.- B
The operator B as a
VectorArray
fromA.source
.- C
The operator C as a
VectorArray
fromA.source
.- R
The operator R as a 2D
NumPy array
orNone
.- S
The operator S as a
VectorArray
fromA.source
orNone
.- trans
Whether the first
Operator
in 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,
VectorArray
fromA.source
.
-
pymor.bindings.slycot.
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
slycot.sb02md
(if E and S areNone
),slycot.sb02od
(if E isNone
and S is notNone
) andslycot.sg03ad
(if E is notNone
), which are dense solvers. Therefore, we assume allOperators
andVectorArrays
can be converted toNumPy arrays
usingto_matrix
andto_numpy
.Parameters
- A
The non-parametric
Operator
A.- E
The non-parametric
Operator
E orNone
.- B
The operator B as a
VectorArray
fromA.source
.- C
The operator C as a
VectorArray
fromA.source
.- R
The operator R as a 2D
NumPy array
orNone
.- S
The operator S as a
VectorArray
fromA.source
orNone
.- 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
fromA.source
.