Skip to content
Merged
Show file tree
Hide file tree
Changes from 1 commit
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
Prev Previous commit
Next Next commit
update docstrings for wrappers of changed methods in SLICOT v5.8
  • Loading branch information
bnavigator committed May 29, 2022
commit 0e83e44368fe14b0e13298ea5f0f55e70159299e
2 changes: 1 addition & 1 deletion slycot/math.py
Original file line number Diff line number Diff line change
Expand Up @@ -29,7 +29,7 @@ def mb03rd(n, A, X=None, jobx='U', sort='N', pmax=1.0, tol=0.0):
To reduce a matrix `A` in real Schur form to a block-diagonal form
using well-conditioned non-orthogonal similarity transformations.
The condition numbers of the transformations used for reduction
are roughly bounded by `pmax`*`pmax`, where `pmax` is a given value.
are roughly bounded by `pmax`, where `pmax` is a given value.
The transformations are optionally postmultiplied in a given
matrix `X`. The real Schur form is optionally ordered, so that
clustered eigenvalues are grouped in the same block.
Expand Down
25 changes: 11 additions & 14 deletions slycot/synthesis.py
Original file line number Diff line number Diff line change
Expand Up @@ -839,7 +839,7 @@ def sb03od(n,m,A,Q,B,dico,fact='N',trans='N',ldwork=None):
set less than or equal to 1 to avoid overflow in X. If matrix B
has full rank then the solution matrix X will be positive-definite
and hence the Cholesky factor U will be nonsingular, but if B is
rank deficient then X may be only positive semi-definite and U
rank deficient, then X may be only positive semi-definite and U
will be singular.

In the case of equation (1) the matrix A must be stable (that
Expand All @@ -850,28 +850,28 @@ def sb03od(n,m,A,Q,B,dico,fact='N',trans='N',ldwork=None):
Parameters
----------
n : int
The order of the matrix A and the number of columns in
matrix op(B). n >= 0.
The order of the matrix A and the number of columns of
the matrix op(B). n >= 0.
m : int
The number of rows in matrix op(B). m >= 0.
A : (n, n) array_like
On entry, the leading n-by-n part of this array must
contain the matrix A. If fact = 'F', then A contains
an upper quasi-triangular matrix S in Schur canonical
form; the elements below the upper Hessenberg part of the
array A are not referenced.
array A are then not referenced.
On exit, the leading n-by-n upper Hessenberg part of this
array contains the upper quasi-triangular matrix S in
Schur canonical form from the Shur factorization of A.
The contents of array A is not modified if fact = 'F'.
The contents of the array A is not modified if fact = 'F'.
Q : (n, n) array_like
On entry, if fact = 'F', then the leading n-by-n part of
this array must contain the orthogonal matrix Q of the
Schur factorization of A.
Otherwise, Q need not be set on entry.
On exit, the leading n-by-n part of this array contains
the orthogonal matrix Q of the Schur factorization of A.
The contents of array Q is not modified if fact = 'F'.
The contents of the array Q is not modified if fact = 'F'.
B : (m, n) array_like
On entry, if trans = 'N', the leading m-by-n part of this
array must contain the coefficient matrix B of the
Expand Down Expand Up @@ -915,7 +915,7 @@ def sb03od(n,m,A,Q,B,dico,fact='N',trans='N',ldwork=None):
The scale factor, scale, set less than or equal to 1 to
prevent the solution overflowing.
w : (n, ) complex ndarray
If fact = 'N', this array contains the eigenvalues of A.
The eigenvalues of A.

Raises
------
Expand All @@ -932,7 +932,7 @@ def sb03od(n,m,A,Q,B,dico,fact='N',trans='N',ldwork=None):
:info = 4:
FACT = 'F' and the Schur factor S supplied in
the array A has two or more consecutive non-zero
elements on the first sub-diagonal, so that there is
elements on the first subdiagonal, so that there is
a block larger than 2-by-2 on the diagonal
:info = 5:
FACT = 'F' and the Schur factor S supplied in
Expand Down Expand Up @@ -2333,9 +2333,6 @@ def sg03bd(n,m,A,E,Q,Z,B,dico,fact='N',trans='N',ldwork=None):
Hessenberg part of the array A are not referenced.
If fact = 'N', then the leading n-by-n part of this
array must contain the matrix A.
On exit, the leading n-by-n part of this array contains
the generalized Schur factor A_s of the matrix A. (A_s is
an upper quasitriangular matrix.)
E : (n, n) array_like
On entry, if fact = 'F', then the leading n-by-n upper
triangular part of this array must contain the
Expand Down Expand Up @@ -2407,9 +2404,9 @@ def sg03bd(n,m,A,E,Q,Z,B,dico,fact='N',trans='N',ldwork=None):
scale : float
The scale factor set to avoid overflow in U.
0 < scale <= 1.
alpha : (n, ) complex ndarray
lambda : (n, ) complex ndarray
If INFO = 0, 3, 5, 6, or 7, then
(alpha(j), j=1,...,n, are the
((j), j=1,...,n, are the
eigenvalues of the matrix pencil A - lambda * E.

Raises
Expand Down Expand Up @@ -2463,7 +2460,7 @@ def sg03bd(n,m,A,E,Q,Z,B,dico,fact='N',trans='N',ldwork=None):
alpha = _np.zeros(n,'complex64')
alpha.real = alphar[0:n]
alpha.imag = alphai[0:n]
return U,scale,alpha/beta
return U, scale, alpha/beta


def sb10fd(n,m,np,ncon,nmeas,gamma,A,B,C,D,tol=0.0,ldwork=None):
Expand Down