pymor.bindings package¶
Submodules¶
fenics module¶
-
class
pymor.bindings.fenics.ComplexifiedFenicsVector(real_part, imag_part)[source]¶ Bases:
pymor.vectorarrays.list.ComplexifiedVectorMethods
amaxaxpy,conj,copy,dofs,inner,norm,norm2,scal,sup_norm,to_numpyAttributes
imag,real
-
class
pymor.bindings.fenics.FenicsMatrixOperator(*args, **kwargs)[source]¶ Bases:
pymor.operators.list.LinearComplexifiedListVectorArrayOperatorBaseWraps a FEniCS matrix as an
Operator.Methods
apply2,as_range_array,as_source_array,as_vector,assemble,d_mu,jacobian,pairwise_apply2,restricted,__matmul__Attributes
parameters,parameters_inherited,parameters_internal,parameters_own,parametric-
_assemble_lincomb(operators, coefficients, identity_shift=0.0, solver_options=None, name=None)[source]¶ Try to assemble a linear combination of the given operators.
Returns a new
Operatorwhich represents the sumc_1*O_1 + ... + c_N*O_N + s*I
where
O_iareOperators,c_i,sscalar coefficients andIthe identity.This method is called in the
assemblemethod ofLincombOperatoron the first of its operators. If an assembly of the given linear combination is possible, e.g. the linear combination of the system matrices of the operators can be formed, then the assembled operator is returned. Otherwise, the method returnsNoneto indicate that assembly is not possible.Parameters
- operators
List of
OperatorsO_iwhose linear combination is formed.- coefficients
List of the corresponding linear coefficients
c_i.- identity_shift
The coefficient
s.- solver_options
solver_optionsfor the assembled operator.- name
Name of the assembled operator.
Returns
The assembled
Operatorif assembly is possible, otherwiseNone.
-
-
class
pymor.bindings.fenics.FenicsOperator(*args, **kwargs)[source]¶ Bases:
pymor.operators.interface.OperatorWraps 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,__matmul__Attributes
parameters,parameters_inherited,parameters_internal,parameters_own,parametric-
apply(U, mu=None)[source]¶ Apply the operator to a
VectorArray.Parameters
- U
VectorArrayof vectors to which the operator is applied.- mu
The
parameter valuesfor which to evaluate the operator.
Returns
VectorArrayof the operator evaluations.
-
jacobian(U, mu=None)[source]¶ Return the operator’s Jacobian as a new
Operator.Parameters
- U
Length 1
VectorArraycontaining the vector for which to compute the Jacobian.- mu
The
parameter valuesfor which to compute the Jacobian.
Returns
Linear
Operatorrepresenting the Jacobian.
-
restricted(dofs)[source]¶ Restrict the operator range to a given set of degrees of freedom.
This method returns a restricted version
restricted_opof the operator along with an arraysource_dofssuch that for anyVectorArrayUinself.sourcethe 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 interpolationwhere 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_dofswill 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 arrayof degrees of freedom in the operatorrangeto which to restrict.
Returns
- restricted_op
The restricted operator as defined above. The operator will have
NumpyVectorSpace(len(source_dofs))assourceandNumpyVectorSpace(len(dofs))asrange.- source_dofs
One-dimensional
NumPy arrayof source degrees of freedom as defined above.
-
-
class
pymor.bindings.fenics.FenicsVector(impl)[source]¶ Bases:
pymor.vectorarrays.list.CopyOnWriteVectorWraps a FEniCS vector to make it usable with ListVectorArray.
Methods
amax,dofs,from_instance,inner,norm,norm2,sup_norm,to_numpyaxpy,copy,scalconjAttributes
imag,real
-
class
pymor.bindings.fenics.FenicsVectorSpace(*args, **kwargs)[source]¶ Bases:
pymor.vectorarrays.list.ComplexifiedListVectorSpaceMethods
real_full_vector,real_make_vector,real_random_vector,real_vector_from_numpy,real_zero_vectorfull_vector,make_vector,random_vector,vector_from_numpy,zero_vectorfrom_numpy,full,make_array,ones,ones_vector,random,space_from_dim,space_from_vector_obj,zerosAttributes
-
complexified_vector_type¶ alias of
ComplexifiedFenicsVector
-
-
class
pymor.bindings.fenics.FenicsVisualizer(*args, **kwargs)[source]¶ Bases:
pymor.core.base.ImmutableObjectVisualize a FEniCS grid function.
Parameters
- space
The
FenicsVectorSpacefor 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).
Methods
Attributes
-
visualize(U, m, title='', legend=None, filename=None, block=True, separate_colorbars=True)[source]¶ Visualize the provided data.
Parameters
- U
VectorArrayof the data to visualize (length must be 1). Alternatively, a tuple ofVectorArrayswhich will be visualized in separate windows. Iffilenameis specified, only oneVectorArraymay 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
Uis a tuple ofVectorArrays,legendhas to be a tuple of the same length.- filename
If specified, write the data to that file.
filenameneeds 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.OperatorMethods
apply2,apply_adjoint,apply_inverse,apply_inverse_adjoint,as_range_array,as_source_array,as_vector,assemble,d_mu,pairwise_apply2,restricted,__matmul__Attributes
parameters,parameters_inherited,parameters_internal,parameters_own,parametric-
apply(U, mu=None)[source]¶ Apply the operator to a
VectorArray.Parameters
- U
VectorArrayof vectors to which the operator is applied.- mu
The
parameter valuesfor which to evaluate the operator.
Returns
VectorArrayof the operator evaluations.
-
jacobian(U, mu=None)[source]¶ Return the operator’s Jacobian as a new
Operator.Parameters
- U
Length 1
VectorArraycontaining the vector for which to compute the Jacobian.- mu
The
parameter valuesfor which to compute the Jacobian.
Returns
Linear
Operatorrepresenting the Jacobian.
-
ngsolve module¶
-
class
pymor.bindings.ngsolve.ComplexifiedNGSolveVector(real_part, imag_part)[source]¶ Bases:
pymor.bindings.ngsolve.NGSolveVectorCommon,pymor.vectorarrays.list.ComplexifiedVectorMethods
amax,dofsaxpy,conj,copy,inner,norm,norm2,scal,sup_norm,to_numpyAttributes
imag,real
-
class
pymor.bindings.ngsolve.NGSolveMatrixOperator(*args, **kwargs)[source]¶ Bases:
pymor.operators.list.LinearComplexifiedListVectorArrayOperatorBaseWraps a NGSolve matrix as an
Operator.Methods
apply2,as_range_array,as_source_array,assemble,d_mu,jacobian,pairwise_apply2,restricted,__matmul__Attributes
parameters,parameters_inherited,parameters_internal,parameters_own,parametric-
_assemble_lincomb(operators, coefficients, identity_shift=0.0, solver_options=None, name=None)[source]¶ Try to assemble a linear combination of the given operators.
Returns a new
Operatorwhich represents the sumc_1*O_1 + ... + c_N*O_N + s*I
where
O_iareOperators,c_i,sscalar coefficients andIthe identity.This method is called in the
assemblemethod ofLincombOperatoron the first of its operators. If an assembly of the given linear combination is possible, e.g. the linear combination of the system matrices of the operators can be formed, then the assembled operator is returned. Otherwise, the method returnsNoneto indicate that assembly is not possible.Parameters
- operators
List of
OperatorsO_iwhose linear combination is formed.- coefficients
List of the corresponding linear coefficients
c_i.- identity_shift
The coefficient
s.- solver_options
solver_optionsfor the assembled operator.- name
Name of the assembled operator.
Returns
The assembled
Operatorif assembly is possible, otherwiseNone.
-
as_vector(copy=True)[source]¶ Return a vector representation of a linear functional or vector operator.
Depending on the operator’s
sourceandrange, this method is equivalent to callingas_range_arrayoras_source_arrayrespectively. The resultingVectorArrayis required to have length 1.Parameters
- mu
The
parameter valuesfor which to return the vector representation.
Returns
- V
VectorArrayof length 1 containing the vector representation.
-
-
class
pymor.bindings.ngsolve.NGSolveVector(impl)[source]¶ Bases:
pymor.bindings.ngsolve.NGSolveVectorCommon,pymor.vectorarrays.list.CopyOnWriteVectorWraps a NGSolve BaseVector to make it usable with ListVectorArray.
Methods
from_instance,inner,norm,norm2,to_numpyamax,dofsaxpy,copy,scalconj,sup_normAttributes
imag,real
-
class
pymor.bindings.ngsolve.NGSolveVectorSpace(*args, **kwargs)[source]¶ Bases:
pymor.vectorarrays.list.ComplexifiedListVectorSpaceMethods
real_make_vector,real_vector_from_numpy,real_zero_vector,space_from_vector_objfull_vector,make_vector,random_vector,real_full_vector,real_random_vector,vector_from_numpy,zero_vectorfrom_numpy,full,make_array,ones,ones_vector,random,space_from_dim,zerosAttributes
complexified_vector_type,dim,value_dim-
complexified_vector_type¶ alias of
ComplexifiedNGSolveVector
-
-
class
pymor.bindings.ngsolve.NGSolveVisualizer(*args, **kwargs)[source]¶ Bases:
pymor.core.base.ImmutableObjectVisualize 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
opto the vectors inVusing PyAMG.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
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
VectorArrayof 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_optionsfor 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.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.
Methods
ainv_apply,apeinv_apply,apex_apply,ax_apply,einv_apply,ex_apply,parameterEquationainv_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
Equationb,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.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.
Methods
ainv_apply,apeinv_apply,apex_apply,ax_apply,einv_apply,ex_apply,parameterEquationainv_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
Equationb,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, 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_optionsandlrnm_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_densefor 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 arrayorNone.- 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_lrcffor 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.- 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.
Defaults
default_solver (see
pymor.core.defaults)
-
pymor.bindings.pymess.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
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 operator R as a 2D
NumPy arrayorNone.- 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, 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 operator 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).- 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.
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
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.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_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 operator A as a 2D
NumPy array.- E
The operator E as a 2D
NumPy arrayorNone.- 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_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.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 operator 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.
-
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 operator 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.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.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_densefor 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 arrayorNone.- 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_lrcffor 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 arraysusingto_matrixand thatB.to_numpyis implemented.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.slycot.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
slycot.sb02md(if E isNone) orslycot.sg03ad(if E is notNone), which are dense solvers. 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 operator 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.
-
pymor.bindings.slycot.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
slycot.sb02md(if E isNone) orslycot.sg03ad(if E is notNone), which are dense solvers. 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 operator 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.