DOUBLE PRECISION routines for triangular (or in some cases quasi-triangular) matrix
dtrcon
USAGE:
rcond, info = NumRu::Lapack.dtrcon( norm, uplo, diag, a, [:usage => usage, :help => help])
FORTRAN MANUAL
SUBROUTINE DTRCON( NORM, UPLO, DIAG, N, A, LDA, RCOND, WORK, IWORK, INFO )
* Purpose
* =======
*
* DTRCON estimates the reciprocal of the condition number of a
* triangular matrix A, in either the 1-norm or the infinity-norm.
*
* The norm of A is computed and an estimate is obtained for
* norm(inv(A)), then the reciprocal of the condition number is
* computed as
* RCOND = 1 / ( norm(A) * norm(inv(A)) ).
*
* Arguments
* =========
*
* NORM (input) CHARACTER*1
* Specifies whether the 1-norm condition number or the
* infinity-norm condition number is required:
* = '1' or 'O': 1-norm;
* = 'I': Infinity-norm.
*
* UPLO (input) CHARACTER*1
* = 'U': A is upper triangular;
* = 'L': A is lower triangular.
*
* DIAG (input) CHARACTER*1
* = 'N': A is non-unit triangular;
* = 'U': A is unit triangular.
*
* N (input) INTEGER
* The order of the matrix A. N >= 0.
*
* A (input) DOUBLE PRECISION array, dimension (LDA,N)
* The triangular matrix A. If UPLO = 'U', the leading N-by-N
* upper triangular part of the array A contains the upper
* triangular matrix, and the strictly lower triangular part of
* A is not referenced. If UPLO = 'L', the leading N-by-N lower
* triangular part of the array A contains the lower triangular
* matrix, and the strictly upper triangular part of A is not
* referenced. If DIAG = 'U', the diagonal elements of A are
* also not referenced and are assumed to be 1.
*
* LDA (input) INTEGER
* The leading dimension of the array A. LDA >= max(1,N).
*
* RCOND (output) DOUBLE PRECISION
* The reciprocal of the condition number of the matrix A,
* computed as RCOND = 1/(norm(A) * norm(inv(A))).
*
* WORK (workspace) DOUBLE PRECISION array, dimension (3*N)
*
* IWORK (workspace) INTEGER array, dimension (N)
*
* INFO (output) INTEGER
* = 0: successful exit
* < 0: if INFO = -i, the i-th argument had an illegal value
*
* =====================================================================
*
go to the page top
dtrevc
USAGE:
m, info, select, vl, vr = NumRu::Lapack.dtrevc( side, howmny, select, t, vl, vr, [:usage => usage, :help => help])
FORTRAN MANUAL
SUBROUTINE DTREVC( SIDE, HOWMNY, SELECT, N, T, LDT, VL, LDVL, VR, LDVR, MM, M, WORK, INFO )
* Purpose
* =======
*
* DTREVC computes some or all of the right and/or left eigenvectors of
* a real upper quasi-triangular matrix T.
* Matrices of this type are produced by the Schur factorization of
* a real general matrix: A = Q*T*Q**T, as computed by DHSEQR.
*
* The right eigenvector x and the left eigenvector y of T corresponding
* to an eigenvalue w are defined by:
*
* T*x = w*x, (y**H)*T = w*(y**H)
*
* where y**H denotes the conjugate transpose of y.
* The eigenvalues are not input to this routine, but are read directly
* from the diagonal blocks of T.
*
* This routine returns the matrices X and/or Y of right and left
* eigenvectors of T, or the products Q*X and/or Q*Y, where Q is an
* input matrix. If Q is the orthogonal factor that reduces a matrix
* A to Schur form T, then Q*X and Q*Y are the matrices of right and
* left eigenvectors of A.
*
* Arguments
* =========
*
* SIDE (input) CHARACTER*1
* = 'R': compute right eigenvectors only;
* = 'L': compute left eigenvectors only;
* = 'B': compute both right and left eigenvectors.
*
* HOWMNY (input) CHARACTER*1
* = 'A': compute all right and/or left eigenvectors;
* = 'B': compute all right and/or left eigenvectors,
* backtransformed by the matrices in VR and/or VL;
* = 'S': compute selected right and/or left eigenvectors,
* as indicated by the logical array SELECT.
*
* SELECT (input/output) LOGICAL array, dimension (N)
* If HOWMNY = 'S', SELECT specifies the eigenvectors to be
* computed.
* If w(j) is a real eigenvalue, the corresponding real
* eigenvector is computed if SELECT(j) is .TRUE..
* If w(j) and w(j+1) are the real and imaginary parts of a
* complex eigenvalue, the corresponding complex eigenvector is
* computed if either SELECT(j) or SELECT(j+1) is .TRUE., and
* on exit SELECT(j) is set to .TRUE. and SELECT(j+1) is set to
* .FALSE..
* Not referenced if HOWMNY = 'A' or 'B'.
*
* N (input) INTEGER
* The order of the matrix T. N >= 0.
*
* T (input) DOUBLE PRECISION array, dimension (LDT,N)
* The upper quasi-triangular matrix T in Schur canonical form.
*
* LDT (input) INTEGER
* The leading dimension of the array T. LDT >= max(1,N).
*
* VL (input/output) DOUBLE PRECISION array, dimension (LDVL,MM)
* On entry, if SIDE = 'L' or 'B' and HOWMNY = 'B', VL must
* contain an N-by-N matrix Q (usually the orthogonal matrix Q
* of Schur vectors returned by DHSEQR).
* On exit, if SIDE = 'L' or 'B', VL contains:
* if HOWMNY = 'A', the matrix Y of left eigenvectors of T;
* if HOWMNY = 'B', the matrix Q*Y;
* if HOWMNY = 'S', the left eigenvectors of T specified by
* SELECT, stored consecutively in the columns
* of VL, in the same order as their
* eigenvalues.
* A complex eigenvector corresponding to a complex eigenvalue
* is stored in two consecutive columns, the first holding the
* real part, and the second the imaginary part.
* Not referenced if SIDE = 'R'.
*
* LDVL (input) INTEGER
* The leading dimension of the array VL. LDVL >= 1, and if
* SIDE = 'L' or 'B', LDVL >= N.
*
* VR (input/output) DOUBLE PRECISION array, dimension (LDVR,MM)
* On entry, if SIDE = 'R' or 'B' and HOWMNY = 'B', VR must
* contain an N-by-N matrix Q (usually the orthogonal matrix Q
* of Schur vectors returned by DHSEQR).
* On exit, if SIDE = 'R' or 'B', VR contains:
* if HOWMNY = 'A', the matrix X of right eigenvectors of T;
* if HOWMNY = 'B', the matrix Q*X;
* if HOWMNY = 'S', the right eigenvectors of T specified by
* SELECT, stored consecutively in the columns
* of VR, in the same order as their
* eigenvalues.
* A complex eigenvector corresponding to a complex eigenvalue
* is stored in two consecutive columns, the first holding the
* real part and the second the imaginary part.
* Not referenced if SIDE = 'L'.
*
* LDVR (input) INTEGER
* The leading dimension of the array VR. LDVR >= 1, and if
* SIDE = 'R' or 'B', LDVR >= N.
*
* MM (input) INTEGER
* The number of columns in the arrays VL and/or VR. MM >= M.
*
* M (output) INTEGER
* The number of columns in the arrays VL and/or VR actually
* used to store the eigenvectors.
* If HOWMNY = 'A' or 'B', M is set to N.
* Each selected real eigenvector occupies one column and each
* selected complex eigenvector occupies two columns.
*
* WORK (workspace) DOUBLE PRECISION array, dimension (3*N)
*
* INFO (output) INTEGER
* = 0: successful exit
* < 0: if INFO = -i, the i-th argument had an illegal value
*
* Further Details
* ===============
*
* The algorithm used in this program is basically backward (forward)
* substitution, with scaling to make the the code robust against
* possible overflow.
*
* Each eigenvector is normalized so that the element of largest
* magnitude has magnitude 1; here the magnitude of a complex number
* (x,y) is taken to be |x| + |y|.
*
* =====================================================================
*
go to the page top
dtrexc
USAGE:
info, t, q, ifst, ilst = NumRu::Lapack.dtrexc( compq, t, q, ifst, ilst, [:usage => usage, :help => help])
FORTRAN MANUAL
SUBROUTINE DTREXC( COMPQ, N, T, LDT, Q, LDQ, IFST, ILST, WORK, INFO )
* Purpose
* =======
*
* DTREXC reorders the real Schur factorization of a real matrix
* A = Q*T*Q**T, so that the diagonal block of T with row index IFST is
* moved to row ILST.
*
* The real Schur form T is reordered by an orthogonal similarity
* transformation Z**T*T*Z, and optionally the matrix Q of Schur vectors
* is updated by postmultiplying it with Z.
*
* T must be in Schur canonical form (as returned by DHSEQR), that is,
* block upper triangular with 1-by-1 and 2-by-2 diagonal blocks; each
* 2-by-2 diagonal block has its diagonal elements equal and its
* off-diagonal elements of opposite sign.
*
* Arguments
* =========
*
* COMPQ (input) CHARACTER*1
* = 'V': update the matrix Q of Schur vectors;
* = 'N': do not update Q.
*
* N (input) INTEGER
* The order of the matrix T. N >= 0.
*
* T (input/output) DOUBLE PRECISION array, dimension (LDT,N)
* On entry, the upper quasi-triangular matrix T, in Schur
* Schur canonical form.
* On exit, the reordered upper quasi-triangular matrix, again
* in Schur canonical form.
*
* LDT (input) INTEGER
* The leading dimension of the array T. LDT >= max(1,N).
*
* Q (input/output) DOUBLE PRECISION array, dimension (LDQ,N)
* On entry, if COMPQ = 'V', the matrix Q of Schur vectors.
* On exit, if COMPQ = 'V', Q has been postmultiplied by the
* orthogonal transformation matrix Z which reorders T.
* If COMPQ = 'N', Q is not referenced.
*
* LDQ (input) INTEGER
* The leading dimension of the array Q. LDQ >= max(1,N).
*
* IFST (input/output) INTEGER
* ILST (input/output) INTEGER
* Specify the reordering of the diagonal blocks of T.
* The block with row index IFST is moved to row ILST, by a
* sequence of transpositions between adjacent blocks.
* On exit, if IFST pointed on entry to the second row of a
* 2-by-2 block, it is changed to point to the first row; ILST
* always points to the first row of the block in its final
* position (which may differ from its input value by +1 or -1).
* 1 <= IFST <= N; 1 <= ILST <= N.
*
* WORK (workspace) DOUBLE PRECISION array, dimension (N)
*
* INFO (output) INTEGER
* = 0: successful exit
* < 0: if INFO = -i, the i-th argument had an illegal value
* = 1: two adjacent blocks were too close to swap (the problem
* is very ill-conditioned); T may have been partially
* reordered, and ILST points to the first row of the
* current position of the block being moved.
*
* =====================================================================
*
go to the page top
dtrrfs
USAGE:
ferr, berr, info = NumRu::Lapack.dtrrfs( uplo, trans, diag, a, b, x, [:usage => usage, :help => help])
FORTRAN MANUAL
SUBROUTINE DTRRFS( UPLO, TRANS, DIAG, N, NRHS, A, LDA, B, LDB, X, LDX, FERR, BERR, WORK, IWORK, INFO )
* Purpose
* =======
*
* DTRRFS provides error bounds and backward error estimates for the
* solution to a system of linear equations with a triangular
* coefficient matrix.
*
* The solution matrix X must be computed by DTRTRS or some other
* means before entering this routine. DTRRFS does not do iterative
* refinement because doing so cannot improve the backward error.
*
* Arguments
* =========
*
* UPLO (input) CHARACTER*1
* = 'U': A is upper triangular;
* = 'L': A is lower triangular.
*
* TRANS (input) CHARACTER*1
* Specifies the form of the system of equations:
* = 'N': A * X = B (No transpose)
* = 'T': A**T * X = B (Transpose)
* = 'C': A**H * X = B (Conjugate transpose = Transpose)
*
* DIAG (input) CHARACTER*1
* = 'N': A is non-unit triangular;
* = 'U': A is unit triangular.
*
* N (input) INTEGER
* The order of the matrix A. N >= 0.
*
* NRHS (input) INTEGER
* The number of right hand sides, i.e., the number of columns
* of the matrices B and X. NRHS >= 0.
*
* A (input) DOUBLE PRECISION array, dimension (LDA,N)
* The triangular matrix A. If UPLO = 'U', the leading N-by-N
* upper triangular part of the array A contains the upper
* triangular matrix, and the strictly lower triangular part of
* A is not referenced. If UPLO = 'L', the leading N-by-N lower
* triangular part of the array A contains the lower triangular
* matrix, and the strictly upper triangular part of A is not
* referenced. If DIAG = 'U', the diagonal elements of A are
* also not referenced and are assumed to be 1.
*
* LDA (input) INTEGER
* The leading dimension of the array A. LDA >= max(1,N).
*
* B (input) DOUBLE PRECISION array, dimension (LDB,NRHS)
* The right hand side matrix B.
*
* LDB (input) INTEGER
* The leading dimension of the array B. LDB >= max(1,N).
*
* X (input) DOUBLE PRECISION array, dimension (LDX,NRHS)
* The solution matrix X.
*
* LDX (input) INTEGER
* The leading dimension of the array X. LDX >= max(1,N).
*
* FERR (output) DOUBLE PRECISION array, dimension (NRHS)
* The estimated forward error bound for each solution vector
* X(j) (the j-th column of the solution matrix X).
* If XTRUE is the true solution corresponding to X(j), FERR(j)
* is an estimated upper bound for the magnitude of the largest
* element in (X(j) - XTRUE) divided by the magnitude of the
* largest element in X(j). The estimate is as reliable as
* the estimate for RCOND, and is almost always a slight
* overestimate of the true error.
*
* BERR (output) DOUBLE PRECISION array, dimension (NRHS)
* The componentwise relative backward error of each solution
* vector X(j) (i.e., the smallest relative change in
* any element of A or B that makes X(j) an exact solution).
*
* WORK (workspace) DOUBLE PRECISION array, dimension (3*N)
*
* IWORK (workspace) INTEGER array, dimension (N)
*
* INFO (output) INTEGER
* = 0: successful exit
* < 0: if INFO = -i, the i-th argument had an illegal value
*
* =====================================================================
*
go to the page top
dtrsen
USAGE:
wr, wi, m, s, sep, work, info, t, q = NumRu::Lapack.dtrsen( job, compq, select, t, q, liwork, [:lwork => lwork, :usage => usage, :help => help])
FORTRAN MANUAL
SUBROUTINE DTRSEN( JOB, COMPQ, SELECT, N, T, LDT, Q, LDQ, WR, WI, M, S, SEP, WORK, LWORK, IWORK, LIWORK, INFO )
* Purpose
* =======
*
* DTRSEN reorders the real Schur factorization of a real matrix
* A = Q*T*Q**T, so that a selected cluster of eigenvalues appears in
* the leading diagonal blocks of the upper quasi-triangular matrix T,
* and the leading columns of Q form an orthonormal basis of the
* corresponding right invariant subspace.
*
* Optionally the routine computes the reciprocal condition numbers of
* the cluster of eigenvalues and/or the invariant subspace.
*
* T must be in Schur canonical form (as returned by DHSEQR), that is,
* block upper triangular with 1-by-1 and 2-by-2 diagonal blocks; each
* 2-by-2 diagonal block has its diagonal elemnts equal and its
* off-diagonal elements of opposite sign.
*
* Arguments
* =========
*
* JOB (input) CHARACTER*1
* Specifies whether condition numbers are required for the
* cluster of eigenvalues (S) or the invariant subspace (SEP):
* = 'N': none;
* = 'E': for eigenvalues only (S);
* = 'V': for invariant subspace only (SEP);
* = 'B': for both eigenvalues and invariant subspace (S and
* SEP).
*
* COMPQ (input) CHARACTER*1
* = 'V': update the matrix Q of Schur vectors;
* = 'N': do not update Q.
*
* SELECT (input) LOGICAL array, dimension (N)
* SELECT specifies the eigenvalues in the selected cluster. To
* select a real eigenvalue w(j), SELECT(j) must be set to
* .TRUE.. To select a complex conjugate pair of eigenvalues
* w(j) and w(j+1), corresponding to a 2-by-2 diagonal block,
* either SELECT(j) or SELECT(j+1) or both must be set to
* .TRUE.; a complex conjugate pair of eigenvalues must be
* either both included in the cluster or both excluded.
*
* N (input) INTEGER
* The order of the matrix T. N >= 0.
*
* T (input/output) DOUBLE PRECISION array, dimension (LDT,N)
* On entry, the upper quasi-triangular matrix T, in Schur
* canonical form.
* On exit, T is overwritten by the reordered matrix T, again in
* Schur canonical form, with the selected eigenvalues in the
* leading diagonal blocks.
*
* LDT (input) INTEGER
* The leading dimension of the array T. LDT >= max(1,N).
*
* Q (input/output) DOUBLE PRECISION array, dimension (LDQ,N)
* On entry, if COMPQ = 'V', the matrix Q of Schur vectors.
* On exit, if COMPQ = 'V', Q has been postmultiplied by the
* orthogonal transformation matrix which reorders T; the
* leading M columns of Q form an orthonormal basis for the
* specified invariant subspace.
* If COMPQ = 'N', Q is not referenced.
*
* LDQ (input) INTEGER
* The leading dimension of the array Q.
* LDQ >= 1; and if COMPQ = 'V', LDQ >= N.
*
* WR (output) DOUBLE PRECISION array, dimension (N)
* WI (output) DOUBLE PRECISION array, dimension (N)
* The real and imaginary parts, respectively, of the reordered
* eigenvalues of T. The eigenvalues are stored in the same
* order as on the diagonal of T, with WR(i) = T(i,i) and, if
* T(i:i+1,i:i+1) is a 2-by-2 diagonal block, WI(i) > 0 and
* WI(i+1) = -WI(i). Note that if a complex eigenvalue is
* sufficiently ill-conditioned, then its value may differ
* significantly from its value before reordering.
*
* M (output) INTEGER
* The dimension of the specified invariant subspace.
* 0 < = M <= N.
*
* S (output) DOUBLE PRECISION
* If JOB = 'E' or 'B', S is a lower bound on the reciprocal
* condition number for the selected cluster of eigenvalues.
* S cannot underestimate the true reciprocal condition number
* by more than a factor of sqrt(N). If M = 0 or N, S = 1.
* If JOB = 'N' or 'V', S is not referenced.
*
* SEP (output) DOUBLE PRECISION
* If JOB = 'V' or 'B', SEP is the estimated reciprocal
* condition number of the specified invariant subspace. If
* M = 0 or N, SEP = norm(T).
* If JOB = 'N' or 'E', SEP is not referenced.
*
* WORK (workspace/output) DOUBLE PRECISION array, dimension (MAX(1,LWORK))
* On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
*
* LWORK (input) INTEGER
* The dimension of the array WORK.
* If JOB = 'N', LWORK >= max(1,N);
* if JOB = 'E', LWORK >= max(1,M*(N-M));
* if JOB = 'V' or 'B', LWORK >= max(1,2*M*(N-M)).
*
* If LWORK = -1, then a workspace query is assumed; the routine
* only calculates the optimal size of the WORK array, returns
* this value as the first entry of the WORK array, and no error
* message related to LWORK is issued by XERBLA.
*
* IWORK (workspace) INTEGER array, dimension (MAX(1,LIWORK))
* On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK.
*
* LIWORK (input) INTEGER
* The dimension of the array IWORK.
* If JOB = 'N' or 'E', LIWORK >= 1;
* if JOB = 'V' or 'B', LIWORK >= max(1,M*(N-M)).
*
* If LIWORK = -1, then a workspace query is assumed; the
* routine only calculates the optimal size of the IWORK array,
* returns this value as the first entry of the IWORK array, and
* no error message related to LIWORK is issued by XERBLA.
*
* INFO (output) INTEGER
* = 0: successful exit
* < 0: if INFO = -i, the i-th argument had an illegal value
* = 1: reordering of T failed because some eigenvalues are too
* close to separate (the problem is very ill-conditioned);
* T may have been partially reordered, and WR and WI
* contain the eigenvalues in the same order as in T; S and
* SEP (if requested) are set to zero.
*
* Further Details
* ===============
*
* DTRSEN first collects the selected eigenvalues by computing an
* orthogonal transformation Z to move them to the top left corner of T.
* In other words, the selected eigenvalues are the eigenvalues of T11
* in:
*
* Z'*T*Z = ( T11 T12 ) n1
* ( 0 T22 ) n2
* n1 n2
*
* where N = n1+n2 and Z' means the transpose of Z. The first n1 columns
* of Z span the specified invariant subspace of T.
*
* If T has been obtained from the real Schur factorization of a matrix
* A = Q*T*Q', then the reordered real Schur factorization of A is given
* by A = (Q*Z)*(Z'*T*Z)*(Q*Z)', and the first n1 columns of Q*Z span
* the corresponding invariant subspace of A.
*
* The reciprocal condition number of the average of the eigenvalues of
* T11 may be returned in S. S lies between 0 (very badly conditioned)
* and 1 (very well conditioned). It is computed as follows. First we
* compute R so that
*
* P = ( I R ) n1
* ( 0 0 ) n2
* n1 n2
*
* is the projector on the invariant subspace associated with T11.
* R is the solution of the Sylvester equation:
*
* T11*R - R*T22 = T12.
*
* Let F-norm(M) denote the Frobenius-norm of M and 2-norm(M) denote
* the two-norm of M. Then S is computed as the lower bound
*
* (1 + F-norm(R)**2)**(-1/2)
*
* on the reciprocal of 2-norm(P), the true reciprocal condition number.
* S cannot underestimate 1 / 2-norm(P) by more than a factor of
* sqrt(N).
*
* An approximate error bound for the computed average of the
* eigenvalues of T11 is
*
* EPS * norm(T) / S
*
* where EPS is the machine precision.
*
* The reciprocal condition number of the right invariant subspace
* spanned by the first n1 columns of Z (or of Q*Z) is returned in SEP.
* SEP is defined as the separation of T11 and T22:
*
* sep( T11, T22 ) = sigma-min( C )
*
* where sigma-min(C) is the smallest singular value of the
* n1*n2-by-n1*n2 matrix
*
* C = kprod( I(n2), T11 ) - kprod( transpose(T22), I(n1) )
*
* I(m) is an m by m identity matrix, and kprod denotes the Kronecker
* product. We estimate sigma-min(C) by the reciprocal of an estimate of
* the 1-norm of inverse(C). The true reciprocal 1-norm of inverse(C)
* cannot differ from sigma-min(C) by more than a factor of sqrt(n1*n2).
*
* When SEP is small, small changes in T can cause large changes in
* the invariant subspace. An approximate bound on the maximum angular
* error in the computed right invariant subspace is
*
* EPS * norm(T) / SEP
*
* =====================================================================
*
go to the page top
dtrsna
USAGE:
s, sep, m, info = NumRu::Lapack.dtrsna( job, howmny, select, t, vl, vr, [:usage => usage, :help => help])
FORTRAN MANUAL
SUBROUTINE DTRSNA( JOB, HOWMNY, SELECT, N, T, LDT, VL, LDVL, VR, LDVR, S, SEP, MM, M, WORK, LDWORK, IWORK, INFO )
* Purpose
* =======
*
* DTRSNA estimates reciprocal condition numbers for specified
* eigenvalues and/or right eigenvectors of a real upper
* quasi-triangular matrix T (or of any matrix Q*T*Q**T with Q
* orthogonal).
*
* T must be in Schur canonical form (as returned by DHSEQR), that is,
* block upper triangular with 1-by-1 and 2-by-2 diagonal blocks; each
* 2-by-2 diagonal block has its diagonal elements equal and its
* off-diagonal elements of opposite sign.
*
* Arguments
* =========
*
* JOB (input) CHARACTER*1
* Specifies whether condition numbers are required for
* eigenvalues (S) or eigenvectors (SEP):
* = 'E': for eigenvalues only (S);
* = 'V': for eigenvectors only (SEP);
* = 'B': for both eigenvalues and eigenvectors (S and SEP).
*
* HOWMNY (input) CHARACTER*1
* = 'A': compute condition numbers for all eigenpairs;
* = 'S': compute condition numbers for selected eigenpairs
* specified by the array SELECT.
*
* SELECT (input) LOGICAL array, dimension (N)
* If HOWMNY = 'S', SELECT specifies the eigenpairs for which
* condition numbers are required. To select condition numbers
* for the eigenpair corresponding to a real eigenvalue w(j),
* SELECT(j) must be set to .TRUE.. To select condition numbers
* corresponding to a complex conjugate pair of eigenvalues w(j)
* and w(j+1), either SELECT(j) or SELECT(j+1) or both, must be
* set to .TRUE..
* If HOWMNY = 'A', SELECT is not referenced.
*
* N (input) INTEGER
* The order of the matrix T. N >= 0.
*
* T (input) DOUBLE PRECISION array, dimension (LDT,N)
* The upper quasi-triangular matrix T, in Schur canonical form.
*
* LDT (input) INTEGER
* The leading dimension of the array T. LDT >= max(1,N).
*
* VL (input) DOUBLE PRECISION array, dimension (LDVL,M)
* If JOB = 'E' or 'B', VL must contain left eigenvectors of T
* (or of any Q*T*Q**T with Q orthogonal), corresponding to the
* eigenpairs specified by HOWMNY and SELECT. The eigenvectors
* must be stored in consecutive columns of VL, as returned by
* DHSEIN or DTREVC.
* If JOB = 'V', VL is not referenced.
*
* LDVL (input) INTEGER
* The leading dimension of the array VL.
* LDVL >= 1; and if JOB = 'E' or 'B', LDVL >= N.
*
* VR (input) DOUBLE PRECISION array, dimension (LDVR,M)
* If JOB = 'E' or 'B', VR must contain right eigenvectors of T
* (or of any Q*T*Q**T with Q orthogonal), corresponding to the
* eigenpairs specified by HOWMNY and SELECT. The eigenvectors
* must be stored in consecutive columns of VR, as returned by
* DHSEIN or DTREVC.
* If JOB = 'V', VR is not referenced.
*
* LDVR (input) INTEGER
* The leading dimension of the array VR.
* LDVR >= 1; and if JOB = 'E' or 'B', LDVR >= N.
*
* S (output) DOUBLE PRECISION array, dimension (MM)
* If JOB = 'E' or 'B', the reciprocal condition numbers of the
* selected eigenvalues, stored in consecutive elements of the
* array. For a complex conjugate pair of eigenvalues two
* consecutive elements of S are set to the same value. Thus
* S(j), SEP(j), and the j-th columns of VL and VR all
* correspond to the same eigenpair (but not in general the
* j-th eigenpair, unless all eigenpairs are selected).
* If JOB = 'V', S is not referenced.
*
* SEP (output) DOUBLE PRECISION array, dimension (MM)
* If JOB = 'V' or 'B', the estimated reciprocal condition
* numbers of the selected eigenvectors, stored in consecutive
* elements of the array. For a complex eigenvector two
* consecutive elements of SEP are set to the same value. If
* the eigenvalues cannot be reordered to compute SEP(j), SEP(j)
* is set to 0; this can only occur when the true value would be
* very small anyway.
* If JOB = 'E', SEP is not referenced.
*
* MM (input) INTEGER
* The number of elements in the arrays S (if JOB = 'E' or 'B')
* and/or SEP (if JOB = 'V' or 'B'). MM >= M.
*
* M (output) INTEGER
* The number of elements of the arrays S and/or SEP actually
* used to store the estimated condition numbers.
* If HOWMNY = 'A', M is set to N.
*
* WORK (workspace) DOUBLE PRECISION array, dimension (LDWORK,N+6)
* If JOB = 'E', WORK is not referenced.
*
* LDWORK (input) INTEGER
* The leading dimension of the array WORK.
* LDWORK >= 1; and if JOB = 'V' or 'B', LDWORK >= N.
*
* IWORK (workspace) INTEGER array, dimension (2*(N-1))
* If JOB = 'E', IWORK is not referenced.
*
* INFO (output) INTEGER
* = 0: successful exit
* < 0: if INFO = -i, the i-th argument had an illegal value
*
* Further Details
* ===============
*
* The reciprocal of the condition number of an eigenvalue lambda is
* defined as
*
* S(lambda) = |v'*u| / (norm(u)*norm(v))
*
* where u and v are the right and left eigenvectors of T corresponding
* to lambda; v' denotes the conjugate-transpose of v, and norm(u)
* denotes the Euclidean norm. These reciprocal condition numbers always
* lie between zero (very badly conditioned) and one (very well
* conditioned). If n = 1, S(lambda) is defined to be 1.
*
* An approximate error bound for a computed eigenvalue W(i) is given by
*
* EPS * norm(T) / S(i)
*
* where EPS is the machine precision.
*
* The reciprocal of the condition number of the right eigenvector u
* corresponding to lambda is defined as follows. Suppose
*
* T = ( lambda c )
* ( 0 T22 )
*
* Then the reciprocal condition number is
*
* SEP( lambda, T22 ) = sigma-min( T22 - lambda*I )
*
* where sigma-min denotes the smallest singular value. We approximate
* the smallest singular value by the reciprocal of an estimate of the
* one-norm of the inverse of T22 - lambda*I. If n = 1, SEP(1) is
* defined to be abs(T(1,1)).
*
* An approximate error bound for a computed right eigenvector VR(i)
* is given by
*
* EPS * norm(T) / SEP(i)
*
* =====================================================================
*
go to the page top
dtrsyl
USAGE:
scale, info, c = NumRu::Lapack.dtrsyl( trana, tranb, isgn, a, b, c, [:usage => usage, :help => help])
FORTRAN MANUAL
SUBROUTINE DTRSYL( TRANA, TRANB, ISGN, M, N, A, LDA, B, LDB, C, LDC, SCALE, INFO )
* Purpose
* =======
*
* DTRSYL solves the real Sylvester matrix equation:
*
* op(A)*X + X*op(B) = scale*C or
* op(A)*X - X*op(B) = scale*C,
*
* where op(A) = A or A**T, and A and B are both upper quasi-
* triangular. A is M-by-M and B is N-by-N; the right hand side C and
* the solution X are M-by-N; and scale is an output scale factor, set
* <= 1 to avoid overflow in X.
*
* A and B must be in Schur canonical form (as returned by DHSEQR), that
* is, block upper triangular with 1-by-1 and 2-by-2 diagonal blocks;
* each 2-by-2 diagonal block has its diagonal elements equal and its
* off-diagonal elements of opposite sign.
*
* Arguments
* =========
*
* TRANA (input) CHARACTER*1
* Specifies the option op(A):
* = 'N': op(A) = A (No transpose)
* = 'T': op(A) = A**T (Transpose)
* = 'C': op(A) = A**H (Conjugate transpose = Transpose)
*
* TRANB (input) CHARACTER*1
* Specifies the option op(B):
* = 'N': op(B) = B (No transpose)
* = 'T': op(B) = B**T (Transpose)
* = 'C': op(B) = B**H (Conjugate transpose = Transpose)
*
* ISGN (input) INTEGER
* Specifies the sign in the equation:
* = +1: solve op(A)*X + X*op(B) = scale*C
* = -1: solve op(A)*X - X*op(B) = scale*C
*
* M (input) INTEGER
* The order of the matrix A, and the number of rows in the
* matrices X and C. M >= 0.
*
* N (input) INTEGER
* The order of the matrix B, and the number of columns in the
* matrices X and C. N >= 0.
*
* A (input) DOUBLE PRECISION array, dimension (LDA,M)
* The upper quasi-triangular matrix A, in Schur canonical form.
*
* LDA (input) INTEGER
* The leading dimension of the array A. LDA >= max(1,M).
*
* B (input) DOUBLE PRECISION array, dimension (LDB,N)
* The upper quasi-triangular matrix B, in Schur canonical form.
*
* LDB (input) INTEGER
* The leading dimension of the array B. LDB >= max(1,N).
*
* C (input/output) DOUBLE PRECISION array, dimension (LDC,N)
* On entry, the M-by-N right hand side matrix C.
* On exit, C is overwritten by the solution matrix X.
*
* LDC (input) INTEGER
* The leading dimension of the array C. LDC >= max(1,M)
*
* SCALE (output) DOUBLE PRECISION
* The scale factor, scale, set <= 1 to avoid overflow in X.
*
* INFO (output) INTEGER
* = 0: successful exit
* < 0: if INFO = -i, the i-th argument had an illegal value
* = 1: A and B have common or very close eigenvalues; perturbed
* values were used to solve the equation (but the matrices
* A and B are unchanged).
*
* =====================================================================
*
go to the page top
dtrti2
USAGE:
info, a = NumRu::Lapack.dtrti2( uplo, diag, a, [:usage => usage, :help => help])
FORTRAN MANUAL
SUBROUTINE DTRTI2( UPLO, DIAG, N, A, LDA, INFO )
* Purpose
* =======
*
* DTRTI2 computes the inverse of a real upper or lower triangular
* matrix.
*
* This is the Level 2 BLAS version of the algorithm.
*
* Arguments
* =========
*
* UPLO (input) CHARACTER*1
* Specifies whether the matrix A is upper or lower triangular.
* = 'U': Upper triangular
* = 'L': Lower triangular
*
* DIAG (input) CHARACTER*1
* Specifies whether or not the matrix A is unit triangular.
* = 'N': Non-unit triangular
* = 'U': Unit triangular
*
* N (input) INTEGER
* The order of the matrix A. N >= 0.
*
* A (input/output) DOUBLE PRECISION array, dimension (LDA,N)
* On entry, the triangular matrix A. If UPLO = 'U', the
* leading n by n upper triangular part of the array A contains
* the upper triangular matrix, and the strictly lower
* triangular part of A is not referenced. If UPLO = 'L', the
* leading n by n lower triangular part of the array A contains
* the lower triangular matrix, and the strictly upper
* triangular part of A is not referenced. If DIAG = 'U', the
* diagonal elements of A are also not referenced and are
* assumed to be 1.
*
* On exit, the (triangular) inverse of the original matrix, in
* the same storage format.
*
* LDA (input) INTEGER
* The leading dimension of the array A. LDA >= max(1,N).
*
* INFO (output) INTEGER
* = 0: successful exit
* < 0: if INFO = -k, the k-th argument had an illegal value
*
* =====================================================================
*
go to the page top
dtrtri
USAGE:
info, a = NumRu::Lapack.dtrtri( uplo, diag, a, [:usage => usage, :help => help])
FORTRAN MANUAL
SUBROUTINE DTRTRI( UPLO, DIAG, N, A, LDA, INFO )
* Purpose
* =======
*
* DTRTRI computes the inverse of a real upper or lower triangular
* matrix A.
*
* This is the Level 3 BLAS version of the algorithm.
*
* Arguments
* =========
*
* UPLO (input) CHARACTER*1
* = 'U': A is upper triangular;
* = 'L': A is lower triangular.
*
* DIAG (input) CHARACTER*1
* = 'N': A is non-unit triangular;
* = 'U': A is unit triangular.
*
* N (input) INTEGER
* The order of the matrix A. N >= 0.
*
* A (input/output) DOUBLE PRECISION array, dimension (LDA,N)
* On entry, the triangular matrix A. If UPLO = 'U', the
* leading N-by-N upper triangular part of the array A contains
* the upper triangular matrix, and the strictly lower
* triangular part of A is not referenced. If UPLO = 'L', the
* leading N-by-N lower triangular part of the array A contains
* the lower triangular matrix, and the strictly upper
* triangular part of A is not referenced. If DIAG = 'U', the
* diagonal elements of A are also not referenced and are
* assumed to be 1.
* On exit, the (triangular) inverse of the original matrix, in
* the same storage format.
*
* LDA (input) INTEGER
* The leading dimension of the array A. LDA >= max(1,N).
*
* INFO (output) INTEGER
* = 0: successful exit
* < 0: if INFO = -i, the i-th argument had an illegal value
* > 0: if INFO = i, A(i,i) is exactly zero. The triangular
* matrix is singular and its inverse can not be computed.
*
* =====================================================================
*
go to the page top
dtrtrs
USAGE:
info, b = NumRu::Lapack.dtrtrs( uplo, trans, diag, a, b, [:usage => usage, :help => help])
FORTRAN MANUAL
SUBROUTINE DTRTRS( UPLO, TRANS, DIAG, N, NRHS, A, LDA, B, LDB, INFO )
* Purpose
* =======
*
* DTRTRS solves a triangular system of the form
*
* A * X = B or A**T * X = B,
*
* where A is a triangular matrix of order N, and B is an N-by-NRHS
* matrix. A check is made to verify that A is nonsingular.
*
* Arguments
* =========
*
* UPLO (input) CHARACTER*1
* = 'U': A is upper triangular;
* = 'L': A is lower triangular.
*
* TRANS (input) CHARACTER*1
* Specifies the form of the system of equations:
* = 'N': A * X = B (No transpose)
* = 'T': A**T * X = B (Transpose)
* = 'C': A**H * X = B (Conjugate transpose = Transpose)
*
* DIAG (input) CHARACTER*1
* = 'N': A is non-unit triangular;
* = 'U': A is unit triangular.
*
* N (input) INTEGER
* The order of the matrix A. N >= 0.
*
* NRHS (input) INTEGER
* The number of right hand sides, i.e., the number of columns
* of the matrix B. NRHS >= 0.
*
* A (input) DOUBLE PRECISION array, dimension (LDA,N)
* The triangular matrix A. If UPLO = 'U', the leading N-by-N
* upper triangular part of the array A contains the upper
* triangular matrix, and the strictly lower triangular part of
* A is not referenced. If UPLO = 'L', the leading N-by-N lower
* triangular part of the array A contains the lower triangular
* matrix, and the strictly upper triangular part of A is not
* referenced. If DIAG = 'U', the diagonal elements of A are
* also not referenced and are assumed to be 1.
*
* LDA (input) INTEGER
* The leading dimension of the array A. LDA >= max(1,N).
*
* B (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS)
* On entry, the right hand side matrix B.
* On exit, if INFO = 0, the solution matrix X.
*
* LDB (input) INTEGER
* The leading dimension of the array B. LDB >= max(1,N).
*
* INFO (output) INTEGER
* = 0: successful exit
* < 0: if INFO = -i, the i-th argument had an illegal value
* > 0: if INFO = i, the i-th diagonal element of A is zero,
* indicating that the matrix is singular and the solutions
* X have not been computed.
*
* =====================================================================
*
go to the page top
dtrttf
USAGE:
arf, info = NumRu::Lapack.dtrttf( transr, uplo, a, [:usage => usage, :help => help])
FORTRAN MANUAL
SUBROUTINE DTRTTF( TRANSR, UPLO, N, A, LDA, ARF, INFO )
* Purpose
* =======
*
* DTRTTF copies a triangular matrix A from standard full format (TR)
* to rectangular full packed format (TF) .
*
* Arguments
* =========
*
* TRANSR (input) CHARACTER*1
* = 'N': ARF in Normal form is wanted;
* = 'T': ARF in Transpose form is wanted.
*
* UPLO (input) CHARACTER*1
* = 'U': Upper triangle of A is stored;
* = 'L': Lower triangle of A is stored.
*
* N (input) INTEGER
* The order of the matrix A. N >= 0.
*
* A (input) DOUBLE PRECISION array, dimension (LDA,N).
* On entry, the triangular matrix A. If UPLO = 'U', the
* leading N-by-N upper triangular part of the array A contains
* the upper triangular matrix, and the strictly lower
* triangular part of A is not referenced. If UPLO = 'L', the
* leading N-by-N lower triangular part of the array A contains
* the lower triangular matrix, and the strictly upper
* triangular part of A is not referenced.
*
* LDA (input) INTEGER
* The leading dimension of the matrix A. LDA >= max(1,N).
*
* ARF (output) DOUBLE PRECISION array, dimension (NT).
* NT=N*(N+1)/2. On exit, the triangular matrix A in RFP format.
*
* INFO (output) INTEGER
* = 0: successful exit
* < 0: if INFO = -i, the i-th argument had an illegal value
*
* Further Details
* ===============
*
* We first consider Rectangular Full Packed (RFP) Format when N is
* even. We give an example where N = 6.
*
* AP is Upper AP is Lower
*
* 00 01 02 03 04 05 00
* 11 12 13 14 15 10 11
* 22 23 24 25 20 21 22
* 33 34 35 30 31 32 33
* 44 45 40 41 42 43 44
* 55 50 51 52 53 54 55
*
*
* Let TRANSR = 'N'. RFP holds AP as follows:
* For UPLO = 'U' the upper trapezoid A(0:5,0:2) consists of the last
* three columns of AP upper. The lower triangle A(4:6,0:2) consists of
* the transpose of the first three columns of AP upper.
* For UPLO = 'L' the lower trapezoid A(1:6,0:2) consists of the first
* three columns of AP lower. The upper triangle A(0:2,0:2) consists of
* the transpose of the last three columns of AP lower.
* This covers the case N even and TRANSR = 'N'.
*
* RFP A RFP A
*
* 03 04 05 33 43 53
* 13 14 15 00 44 54
* 23 24 25 10 11 55
* 33 34 35 20 21 22
* 00 44 45 30 31 32
* 01 11 55 40 41 42
* 02 12 22 50 51 52
*
* Now let TRANSR = 'T'. RFP A in both UPLO cases is just the
* transpose of RFP A above. One therefore gets:
*
*
* RFP A RFP A
*
* 03 13 23 33 00 01 02 33 00 10 20 30 40 50
* 04 14 24 34 44 11 12 43 44 11 21 31 41 51
* 05 15 25 35 45 55 22 53 54 55 22 32 42 52
*
*
* We then consider Rectangular Full Packed (RFP) Format when N is
* odd. We give an example where N = 5.
*
* AP is Upper AP is Lower
*
* 00 01 02 03 04 00
* 11 12 13 14 10 11
* 22 23 24 20 21 22
* 33 34 30 31 32 33
* 44 40 41 42 43 44
*
*
* Let TRANSR = 'N'. RFP holds AP as follows:
* For UPLO = 'U' the upper trapezoid A(0:4,0:2) consists of the last
* three columns of AP upper. The lower triangle A(3:4,0:1) consists of
* the transpose of the first two columns of AP upper.
* For UPLO = 'L' the lower trapezoid A(0:4,0:2) consists of the first
* three columns of AP lower. The upper triangle A(0:1,1:2) consists of
* the transpose of the last two columns of AP lower.
* This covers the case N odd and TRANSR = 'N'.
*
* RFP A RFP A
*
* 02 03 04 00 33 43
* 12 13 14 10 11 44
* 22 23 24 20 21 22
* 00 33 34 30 31 32
* 01 11 44 40 41 42
*
* Now let TRANSR = 'T'. RFP A in both UPLO cases is just the
* transpose of RFP A above. One therefore gets:
*
* RFP A RFP A
*
* 02 12 22 00 01 00 10 20 30 40 50
* 03 13 23 33 11 33 11 21 31 41 51
* 04 14 24 34 44 43 44 22 32 42 52
*
* Reference
* =========
*
* =====================================================================
*
* ..
* .. Local Scalars ..
LOGICAL LOWER, NISODD, NORMALTRANSR
INTEGER I, IJ, J, K, L, N1, N2, NT, NX2, NP1X2
* ..
* .. External Functions ..
LOGICAL LSAME
EXTERNAL LSAME
* ..
* .. External Subroutines ..
EXTERNAL XERBLA
* ..
* .. Intrinsic Functions ..
INTRINSIC MAX, MOD
* ..
go to the page top
dtrttp
USAGE:
ap, info = NumRu::Lapack.dtrttp( uplo, a, [:usage => usage, :help => help])
FORTRAN MANUAL
SUBROUTINE DTRTTP( UPLO, N, A, LDA, AP, INFO )
* Purpose
* =======
*
* DTRTTP copies a triangular matrix A from full format (TR) to standard
* packed format (TP).
*
* Arguments
* =========
*
* UPLO (input) CHARACTER*1
* = 'U': A is upper triangular.
* = 'L': A is lower triangular.
*
* N (input) INTEGER
* The order of the matrices AP and A. N >= 0.
*
* A (input) DOUBLE PRECISION array, dimension (LDA,N)
* On exit, the triangular matrix A. If UPLO = 'U', the leading
* N-by-N upper triangular part of A contains the upper
* triangular part of the matrix A, and the strictly lower
* triangular part of A is not referenced. If UPLO = 'L', the
* leading N-by-N lower triangular part of A contains the lower
* triangular part of the matrix A, and the strictly upper
* triangular part of A is not referenced.
*
* LDA (input) INTEGER
* The leading dimension of the array A. LDA >= max(1,N).
*
* AP (output) DOUBLE PRECISION array, dimension (N*(N+1)/2
* On exit, the upper or lower triangular matrix A, packed
* columnwise in a linear array. The j-th column of A is stored
* in the array AP as follows:
* if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
* if UPLO = 'L', AP(i + (j-1)*(2n-j)/2) = A(i,j) for j<=i<=n.
*
* INFO (output) INTEGER
* = 0: successful exit
* < 0: if INFO = -i, the i-th argument had an illegal value
*
* =====================================================================
*
go to the page top
back to matrix types
back to data types