LAPACK
3.4.2
LAPACK: Linear Algebra PACKage
|
Functions/Subroutines | |
subroutine | zgebak (JOB, SIDE, N, ILO, IHI, SCALE, M, V, LDV, INFO) |
ZGEBAK | |
subroutine | zgebal (JOB, N, A, LDA, ILO, IHI, SCALE, INFO) |
ZGEBAL | |
subroutine | zgebd2 (M, N, A, LDA, D, E, TAUQ, TAUP, WORK, INFO) |
ZGEBD2 reduces a general matrix to bidiagonal form using an unblocked algorithm. | |
subroutine | zgebrd (M, N, A, LDA, D, E, TAUQ, TAUP, WORK, LWORK, INFO) |
ZGEBRD | |
subroutine | zgecon (NORM, N, A, LDA, ANORM, RCOND, WORK, RWORK, INFO) |
ZGECON | |
subroutine | zgeequ (M, N, A, LDA, R, C, ROWCND, COLCND, AMAX, INFO) |
ZGEEQU | |
subroutine | zgeequb (M, N, A, LDA, R, C, ROWCND, COLCND, AMAX, INFO) |
ZGEEQUB | |
subroutine | zgehd2 (N, ILO, IHI, A, LDA, TAU, WORK, INFO) |
ZGEHD2 reduces a general square matrix to upper Hessenberg form using an unblocked algorithm. | |
subroutine | zgehrd (N, ILO, IHI, A, LDA, TAU, WORK, LWORK, INFO) |
ZGEHRD | |
subroutine | zgelq2 (M, N, A, LDA, TAU, WORK, INFO) |
ZGELQ2 computes the LQ factorization of a general rectangular matrix using an unblocked algorithm. | |
subroutine | zgelqf (M, N, A, LDA, TAU, WORK, LWORK, INFO) |
ZGELQF | |
subroutine | zgemqrt (SIDE, TRANS, M, N, K, NB, V, LDV, T, LDT, C, LDC, WORK, INFO) |
ZGEMQRT | |
subroutine | zgeql2 (M, N, A, LDA, TAU, WORK, INFO) |
ZGEQL2 computes the QL factorization of a general rectangular matrix using an unblocked algorithm. | |
subroutine | zgeqlf (M, N, A, LDA, TAU, WORK, LWORK, INFO) |
ZGEQLF | |
subroutine | zgeqp3 (M, N, A, LDA, JPVT, TAU, WORK, LWORK, RWORK, INFO) |
ZGEQP3 | |
subroutine | zgeqpf (M, N, A, LDA, JPVT, TAU, WORK, RWORK, INFO) |
ZGEQPF | |
subroutine | zgeqr2 (M, N, A, LDA, TAU, WORK, INFO) |
ZGEQR2 computes the QR factorization of a general rectangular matrix using an unblocked algorithm. | |
subroutine | zgeqr2p (M, N, A, LDA, TAU, WORK, INFO) |
ZGEQR2P computes the QR factorization of a general rectangular matrix with non-negative diagonal elements using an unblocked algorithm. | |
subroutine | zgeqrf (M, N, A, LDA, TAU, WORK, LWORK, INFO) |
ZGEQRF | |
subroutine | zgeqrfp (M, N, A, LDA, TAU, WORK, LWORK, INFO) |
ZGEQRFP | |
subroutine | zgeqrt (M, N, NB, A, LDA, T, LDT, WORK, INFO) |
ZGEQRT | |
subroutine | zgeqrt2 (M, N, A, LDA, T, LDT, INFO) |
ZGEQRT2 computes a QR factorization of a general real or complex matrix using the compact WY representation of Q. | |
recursive subroutine | zgeqrt3 (M, N, A, LDA, T, LDT, INFO) |
ZGEQRT3 recursively computes a QR factorization of a general real or complex matrix using the compact WY representation of Q. | |
subroutine | zgerfs (TRANS, N, NRHS, A, LDA, AF, LDAF, IPIV, B, LDB, X, LDX, FERR, BERR, WORK, RWORK, INFO) |
ZGERFS | |
subroutine | zgerfsx (TRANS, EQUED, N, NRHS, A, LDA, AF, LDAF, IPIV, R, C, B, LDB, X, LDX, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP, NPARAMS, PARAMS, WORK, RWORK, INFO) |
ZGERFSX | |
subroutine | zgerq2 (M, N, A, LDA, TAU, WORK, INFO) |
ZGERQ2 computes the RQ factorization of a general rectangular matrix using an unblocked algorithm. | |
subroutine | zgerqf (M, N, A, LDA, TAU, WORK, LWORK, INFO) |
ZGERQF | |
subroutine | zgetf2 (M, N, A, LDA, IPIV, INFO) |
ZGETF2 computes the LU factorization of a general m-by-n matrix using partial pivoting with row interchanges (unblocked algorithm). | |
subroutine | zgetrf (M, N, A, LDA, IPIV, INFO) |
ZGETRF | |
subroutine | zgetri (N, A, LDA, IPIV, WORK, LWORK, INFO) |
ZGETRI | |
subroutine | zgetrs (TRANS, N, NRHS, A, LDA, IPIV, B, LDB, INFO) |
ZGETRS | |
subroutine | zhgeqz (JOB, COMPQ, COMPZ, N, ILO, IHI, H, LDH, T, LDT, ALPHA, BETA, Q, LDQ, Z, LDZ, WORK, LWORK, RWORK, INFO) |
ZHGEQZ | |
subroutine | zla_geamv (TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY) |
ZLA_GEAMV computes a matrix-vector product using a general matrix to calculate error bounds. | |
DOUBLE PRECISION function | zla_gercond_c (TRANS, N, A, LDA, AF, LDAF, IPIV, C, CAPPLY, INFO, WORK, RWORK) |
ZLA_GERCOND_C computes the infinity norm condition number of op(A)*inv(diag(c)) for general matrices. | |
DOUBLE PRECISION function | zla_gercond_x (TRANS, N, A, LDA, AF, LDAF, IPIV, X, INFO, WORK, RWORK) |
ZLA_GERCOND_X computes the infinity norm condition number of op(A)*diag(x) for general matrices. | |
subroutine | zla_gerfsx_extended (PREC_TYPE, TRANS_TYPE, N, NRHS, A, LDA, AF, LDAF, IPIV, COLEQU, C, B, LDB, Y, LDY, BERR_OUT, N_NORMS, ERRS_N, ERRS_C, RES, AYB, DY, Y_TAIL, RCOND, ITHRESH, RTHRESH, DZ_UB, IGNORE_CWISE, INFO) |
ZLA_GERFSX_EXTENDED | |
DOUBLE PRECISION function | zla_gerpvgrw (N, NCOLS, A, LDA, AF, LDAF) |
ZLA_GERPVGRW multiplies a square real matrix by a complex matrix. | |
subroutine | ztgevc (SIDE, HOWMNY, SELECT, N, S, LDS, P, LDP, VL, LDVL, VR, LDVR, MM, M, WORK, RWORK, INFO) |
ZTGEVC | |
subroutine | ztgexc (WANTQ, WANTZ, N, A, LDA, B, LDB, Q, LDQ, Z, LDZ, IFST, ILST, INFO) |
ZTGEXC |
This is the group of complex16 computational functions for GE matrices
subroutine zgebak | ( | character | JOB, |
character | SIDE, | ||
integer | N, | ||
integer | ILO, | ||
integer | IHI, | ||
double precision, dimension( * ) | SCALE, | ||
integer | M, | ||
complex*16, dimension( ldv, * ) | V, | ||
integer | LDV, | ||
integer | INFO | ||
) |
ZGEBAK
Download ZGEBAK + dependencies [TGZ] [ZIP] [TXT]ZGEBAK forms the right or left eigenvectors of a complex general matrix by backward transformation on the computed eigenvectors of the balanced matrix output by ZGEBAL.
[in] | JOB | JOB is CHARACTER*1 Specifies the type of backward transformation required: = 'N', do nothing, return immediately; = 'P', do backward transformation for permutation only; = 'S', do backward transformation for scaling only; = 'B', do backward transformations for both permutation and scaling. JOB must be the same as the argument JOB supplied to ZGEBAL. |
[in] | SIDE | SIDE is CHARACTER*1 = 'R': V contains right eigenvectors; = 'L': V contains left eigenvectors. |
[in] | N | N is INTEGER The number of rows of the matrix V. N >= 0. |
[in] | ILO | ILO is INTEGER |
[in] | IHI | IHI is INTEGER The integers ILO and IHI determined by ZGEBAL. 1 <= ILO <= IHI <= N, if N > 0; ILO=1 and IHI=0, if N=0. |
[in] | SCALE | SCALE is DOUBLE PRECISION array, dimension (N) Details of the permutation and scaling factors, as returned by ZGEBAL. |
[in] | M | M is INTEGER The number of columns of the matrix V. M >= 0. |
[in,out] | V | V is COMPLEX*16 array, dimension (LDV,M) On entry, the matrix of right or left eigenvectors to be transformed, as returned by ZHSEIN or ZTREVC. On exit, V is overwritten by the transformed eigenvectors. |
[in] | LDV | LDV is INTEGER The leading dimension of the array V. LDV >= max(1,N). |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value. |
Definition at line 131 of file zgebak.f.
subroutine zgebal | ( | character | JOB, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
integer | ILO, | ||
integer | IHI, | ||
double precision, dimension( * ) | SCALE, | ||
integer | INFO | ||
) |
ZGEBAL
Download ZGEBAL + dependencies [TGZ] [ZIP] [TXT]ZGEBAL balances a general complex matrix A. This involves, first, permuting A by a similarity transformation to isolate eigenvalues in the first 1 to ILO-1 and last IHI+1 to N elements on the diagonal; and second, applying a diagonal similarity transformation to rows and columns ILO to IHI to make the rows and columns as close in norm as possible. Both steps are optional. Balancing may reduce the 1-norm of the matrix, and improve the accuracy of the computed eigenvalues and/or eigenvectors.
[in] | JOB | JOB is CHARACTER*1 Specifies the operations to be performed on A: = 'N': none: simply set ILO = 1, IHI = N, SCALE(I) = 1.0 for i = 1,...,N; = 'P': permute only; = 'S': scale only; = 'B': both permute and scale. |
[in] | N | N is INTEGER The order of the matrix A. N >= 0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the input matrix A. On exit, A is overwritten by the balanced matrix. If JOB = 'N', A is not referenced. See Further Details. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[out] | ILO | |
[out] | IHI | ILO and IHI are set to INTEGER such that on exit A(i,j) = 0 if i > j and j = 1,...,ILO-1 or I = IHI+1,...,N. If JOB = 'N' or 'S', ILO = 1 and IHI = N. |
[out] | SCALE | SCALE is DOUBLE PRECISION array, dimension (N) Details of the permutations and scaling factors applied to A. If P(j) is the index of the row and column interchanged with row and column j and D(j) is the scaling factor applied to row and column j, then SCALE(j) = P(j) for j = 1,...,ILO-1 = D(j) for j = ILO,...,IHI = P(j) for j = IHI+1,...,N. The order in which the interchanges are made is N to IHI+1, then 1 to ILO-1. |
[out] | INFO | INFO is INTEGER = 0: successful exit. < 0: if INFO = -i, the i-th argument had an illegal value. |
The permutations consist of row and column interchanges which put the matrix in the form ( T1 X Y ) P A P = ( 0 B Z ) ( 0 0 T2 ) where T1 and T2 are upper triangular matrices whose eigenvalues lie along the diagonal. The column indices ILO and IHI mark the starting and ending columns of the submatrix B. Balancing consists of applying a diagonal similarity transformation inv(D) * B * D to make the 1-norms of each row of B and its corresponding column nearly equal. The output matrix is ( T1 X*D Y ) ( 0 inv(D)*B*D inv(D)*Z ). ( 0 0 T2 ) Information about the permutations P and the diagonal matrix D is returned in the vector SCALE. This subroutine is based on the EISPACK routine CBAL. Modified by Tzu-Yi Chen, Computer Science Division, University of California at Berkeley, USA
Definition at line 161 of file zgebal.f.
subroutine zgebd2 | ( | integer | M, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
double precision, dimension( * ) | D, | ||
double precision, dimension( * ) | E, | ||
complex*16, dimension( * ) | TAUQ, | ||
complex*16, dimension( * ) | TAUP, | ||
complex*16, dimension( * ) | WORK, | ||
integer | INFO | ||
) |
ZGEBD2 reduces a general matrix to bidiagonal form using an unblocked algorithm.
Download ZGEBD2 + dependencies [TGZ] [ZIP] [TXT]ZGEBD2 reduces a complex general m by n matrix A to upper or lower real bidiagonal form B by a unitary transformation: Q**H * A * P = B. If m >= n, B is upper bidiagonal; if m < n, B is lower bidiagonal.
[in] | M | M is INTEGER The number of rows in the matrix A. M >= 0. |
[in] | N | N is INTEGER The number of columns in the matrix A. N >= 0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the m by n general matrix to be reduced. On exit, if m >= n, the diagonal and the first superdiagonal are overwritten with the upper bidiagonal matrix B; the elements below the diagonal, with the array TAUQ, represent the unitary matrix Q as a product of elementary reflectors, and the elements above the first superdiagonal, with the array TAUP, represent the unitary matrix P as a product of elementary reflectors; if m < n, the diagonal and the first subdiagonal are overwritten with the lower bidiagonal matrix B; the elements below the first subdiagonal, with the array TAUQ, represent the unitary matrix Q as a product of elementary reflectors, and the elements above the diagonal, with the array TAUP, represent the unitary matrix P as a product of elementary reflectors. See Further Details. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | D | D is DOUBLE PRECISION array, dimension (min(M,N)) The diagonal elements of the bidiagonal matrix B: D(i) = A(i,i). |
[out] | E | E is DOUBLE PRECISION array, dimension (min(M,N)-1) The off-diagonal elements of the bidiagonal matrix B: if m >= n, E(i) = A(i,i+1) for i = 1,2,...,n-1; if m < n, E(i) = A(i+1,i) for i = 1,2,...,m-1. |
[out] | TAUQ | TAUQ is COMPLEX*16 array dimension (min(M,N)) The scalar factors of the elementary reflectors which represent the unitary matrix Q. See Further Details. |
[out] | TAUP | TAUP is COMPLEX*16 array, dimension (min(M,N)) The scalar factors of the elementary reflectors which represent the unitary matrix P. See Further Details. |
[out] | WORK | WORK is COMPLEX*16 array, dimension (max(M,N)) |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value. |
The matrices Q and P are represented as products of elementary reflectors: If m >= n, Q = H(1) H(2) . . . H(n) and P = G(1) G(2) . . . G(n-1) Each H(i) and G(i) has the form: H(i) = I - tauq * v * v**H and G(i) = I - taup * u * u**H where tauq and taup are complex scalars, and v and u are complex vectors; v(1:i-1) = 0, v(i) = 1, and v(i+1:m) is stored on exit in A(i+1:m,i); u(1:i) = 0, u(i+1) = 1, and u(i+2:n) is stored on exit in A(i,i+2:n); tauq is stored in TAUQ(i) and taup in TAUP(i). If m < n, Q = H(1) H(2) . . . H(m-1) and P = G(1) G(2) . . . G(m) Each H(i) and G(i) has the form: H(i) = I - tauq * v * v**H and G(i) = I - taup * u * u**H where tauq and taup are complex scalars, v and u are complex vectors; v(1:i) = 0, v(i+1) = 1, and v(i+2:m) is stored on exit in A(i+2:m,i); u(1:i-1) = 0, u(i) = 1, and u(i+1:n) is stored on exit in A(i,i+1:n); tauq is stored in TAUQ(i) and taup in TAUP(i). The contents of A on exit are illustrated by the following examples: m = 6 and n = 5 (m > n): m = 5 and n = 6 (m < n): ( d e u1 u1 u1 ) ( d u1 u1 u1 u1 u1 ) ( v1 d e u2 u2 ) ( e d u2 u2 u2 u2 ) ( v1 v2 d e u3 ) ( v1 e d u3 u3 u3 ) ( v1 v2 v3 d e ) ( v1 v2 e d u4 u4 ) ( v1 v2 v3 v4 d ) ( v1 v2 v3 e d u5 ) ( v1 v2 v3 v4 v5 ) where d and e denote diagonal and off-diagonal elements of B, vi denotes an element of the vector defining H(i), and ui an element of the vector defining G(i).
Definition at line 190 of file zgebd2.f.
subroutine zgebrd | ( | integer | M, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
double precision, dimension( * ) | D, | ||
double precision, dimension( * ) | E, | ||
complex*16, dimension( * ) | TAUQ, | ||
complex*16, dimension( * ) | TAUP, | ||
complex*16, dimension( * ) | WORK, | ||
integer | LWORK, | ||
integer | INFO | ||
) |
ZGEBRD
Download ZGEBRD + dependencies [TGZ] [ZIP] [TXT]ZGEBRD reduces a general complex M-by-N matrix A to upper or lower bidiagonal form B by a unitary transformation: Q**H * A * P = B. If m >= n, B is upper bidiagonal; if m < n, B is lower bidiagonal.
[in] | M | M is INTEGER The number of rows in the matrix A. M >= 0. |
[in] | N | N is INTEGER The number of columns in the matrix A. N >= 0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the M-by-N general matrix to be reduced. On exit, if m >= n, the diagonal and the first superdiagonal are overwritten with the upper bidiagonal matrix B; the elements below the diagonal, with the array TAUQ, represent the unitary matrix Q as a product of elementary reflectors, and the elements above the first superdiagonal, with the array TAUP, represent the unitary matrix P as a product of elementary reflectors; if m < n, the diagonal and the first subdiagonal are overwritten with the lower bidiagonal matrix B; the elements below the first subdiagonal, with the array TAUQ, represent the unitary matrix Q as a product of elementary reflectors, and the elements above the diagonal, with the array TAUP, represent the unitary matrix P as a product of elementary reflectors. See Further Details. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | D | D is DOUBLE PRECISION array, dimension (min(M,N)) The diagonal elements of the bidiagonal matrix B: D(i) = A(i,i). |
[out] | E | E is DOUBLE PRECISION array, dimension (min(M,N)-1) The off-diagonal elements of the bidiagonal matrix B: if m >= n, E(i) = A(i,i+1) for i = 1,2,...,n-1; if m < n, E(i) = A(i+1,i) for i = 1,2,...,m-1. |
[out] | TAUQ | TAUQ is COMPLEX*16 array dimension (min(M,N)) The scalar factors of the elementary reflectors which represent the unitary matrix Q. See Further Details. |
[out] | TAUP | TAUP is COMPLEX*16 array, dimension (min(M,N)) The scalar factors of the elementary reflectors which represent the unitary matrix P. See Further Details. |
[out] | WORK | WORK is COMPLEX*16 array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK(1) returns the optimal LWORK. |
[in] | LWORK | LWORK is INTEGER The length of the array WORK. LWORK >= max(1,M,N). For optimum performance LWORK >= (M+N)*NB, where NB is the optimal blocksize. 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. |
[out] | INFO | INFO is INTEGER = 0: successful exit. < 0: if INFO = -i, the i-th argument had an illegal value. |
The matrices Q and P are represented as products of elementary reflectors: If m >= n, Q = H(1) H(2) . . . H(n) and P = G(1) G(2) . . . G(n-1) Each H(i) and G(i) has the form: H(i) = I - tauq * v * v**H and G(i) = I - taup * u * u**H where tauq and taup are complex scalars, and v and u are complex vectors; v(1:i-1) = 0, v(i) = 1, and v(i+1:m) is stored on exit in A(i+1:m,i); u(1:i) = 0, u(i+1) = 1, and u(i+2:n) is stored on exit in A(i,i+2:n); tauq is stored in TAUQ(i) and taup in TAUP(i). If m < n, Q = H(1) H(2) . . . H(m-1) and P = G(1) G(2) . . . G(m) Each H(i) and G(i) has the form: H(i) = I - tauq * v * v**H and G(i) = I - taup * u * u**H where tauq and taup are complex scalars, and v and u are complex vectors; v(1:i) = 0, v(i+1) = 1, and v(i+2:m) is stored on exit in A(i+2:m,i); u(1:i-1) = 0, u(i) = 1, and u(i+1:n) is stored on exit in A(i,i+1:n); tauq is stored in TAUQ(i) and taup in TAUP(i). The contents of A on exit are illustrated by the following examples: m = 6 and n = 5 (m > n): m = 5 and n = 6 (m < n): ( d e u1 u1 u1 ) ( d u1 u1 u1 u1 u1 ) ( v1 d e u2 u2 ) ( e d u2 u2 u2 u2 ) ( v1 v2 d e u3 ) ( v1 e d u3 u3 u3 ) ( v1 v2 v3 d e ) ( v1 v2 e d u4 u4 ) ( v1 v2 v3 v4 d ) ( v1 v2 v3 e d u5 ) ( v1 v2 v3 v4 v5 ) where d and e denote diagonal and off-diagonal elements of B, vi denotes an element of the vector defining H(i), and ui an element of the vector defining G(i).
Definition at line 205 of file zgebrd.f.
subroutine zgecon | ( | character | NORM, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
double precision | ANORM, | ||
double precision | RCOND, | ||
complex*16, dimension( * ) | WORK, | ||
double precision, dimension( * ) | RWORK, | ||
integer | INFO | ||
) |
ZGECON
Download ZGECON + dependencies [TGZ] [ZIP] [TXT]ZGECON estimates the reciprocal of the condition number of a general complex matrix A, in either the 1-norm or the infinity-norm, using the LU factorization computed by ZGETRF. An estimate is obtained for norm(inv(A)), and the reciprocal of the condition number is computed as RCOND = 1 / ( norm(A) * norm(inv(A)) ).
[in] | NORM | NORM is 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. |
[in] | N | N is INTEGER The order of the matrix A. N >= 0. |
[in] | A | A is COMPLEX*16 array, dimension (LDA,N) The factors L and U from the factorization A = P*L*U as computed by ZGETRF. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in] | ANORM | ANORM is DOUBLE PRECISION If NORM = '1' or 'O', the 1-norm of the original matrix A. If NORM = 'I', the infinity-norm of the original matrix A. |
[out] | RCOND | RCOND is DOUBLE PRECISION The reciprocal of the condition number of the matrix A, computed as RCOND = 1/(norm(A) * norm(inv(A))). |
[out] | WORK | WORK is COMPLEX*16 array, dimension (2*N) |
[out] | RWORK | RWORK is DOUBLE PRECISION array, dimension (2*N) |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value |
Definition at line 124 of file zgecon.f.
subroutine zgeequ | ( | integer | M, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
double precision, dimension( * ) | R, | ||
double precision, dimension( * ) | C, | ||
double precision | ROWCND, | ||
double precision | COLCND, | ||
double precision | AMAX, | ||
integer | INFO | ||
) |
ZGEEQU
Download ZGEEQU + dependencies [TGZ] [ZIP] [TXT]ZGEEQU computes row and column scalings intended to equilibrate an M-by-N matrix A and reduce its condition number. R returns the row scale factors and C the column scale factors, chosen to try to make the largest element in each row and column of the matrix B with elements B(i,j)=R(i)*A(i,j)*C(j) have absolute value 1. R(i) and C(j) are restricted to be between SMLNUM = smallest safe number and BIGNUM = largest safe number. Use of these scaling factors is not guaranteed to reduce the condition number of A but works well in practice.
[in] | M | M is INTEGER The number of rows of the matrix A. M >= 0. |
[in] | N | N is INTEGER The number of columns of the matrix A. N >= 0. |
[in] | A | A is COMPLEX*16 array, dimension (LDA,N) The M-by-N matrix whose equilibration factors are to be computed. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | R | R is DOUBLE PRECISION array, dimension (M) If INFO = 0 or INFO > M, R contains the row scale factors for A. |
[out] | C | C is DOUBLE PRECISION array, dimension (N) If INFO = 0, C contains the column scale factors for A. |
[out] | ROWCND | ROWCND is DOUBLE PRECISION If INFO = 0 or INFO > M, ROWCND contains the ratio of the smallest R(i) to the largest R(i). If ROWCND >= 0.1 and AMAX is neither too large nor too small, it is not worth scaling by R. |
[out] | COLCND | COLCND is DOUBLE PRECISION If INFO = 0, COLCND contains the ratio of the smallest C(i) to the largest C(i). If COLCND >= 0.1, it is not worth scaling by C. |
[out] | AMAX | AMAX is DOUBLE PRECISION Absolute value of largest matrix element. If AMAX is very close to overflow or very close to underflow, the matrix should be scaled. |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value > 0: if INFO = i, and i is <= M: the i-th row of A is exactly zero > M: the (i-M)-th column of A is exactly zero |
Definition at line 140 of file zgeequ.f.
subroutine zgeequb | ( | integer | M, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
double precision, dimension( * ) | R, | ||
double precision, dimension( * ) | C, | ||
double precision | ROWCND, | ||
double precision | COLCND, | ||
double precision | AMAX, | ||
integer | INFO | ||
) |
ZGEEQUB
Download ZGEEQUB + dependencies [TGZ] [ZIP] [TXT]ZGEEQUB computes row and column scalings intended to equilibrate an M-by-N matrix A and reduce its condition number. R returns the row scale factors and C the column scale factors, chosen to try to make the largest element in each row and column of the matrix B with elements B(i,j)=R(i)*A(i,j)*C(j) have an absolute value of at most the radix. R(i) and C(j) are restricted to be a power of the radix between SMLNUM = smallest safe number and BIGNUM = largest safe number. Use of these scaling factors is not guaranteed to reduce the condition number of A but works well in practice. This routine differs from ZGEEQU by restricting the scaling factors to a power of the radix. Baring over- and underflow, scaling by these factors introduces no additional rounding errors. However, the scaled entries' magnitured are no longer approximately 1 but lie between sqrt(radix) and 1/sqrt(radix).
[in] | M | M is INTEGER The number of rows of the matrix A. M >= 0. |
[in] | N | N is INTEGER The number of columns of the matrix A. N >= 0. |
[in] | A | A is COMPLEX*16 array, dimension (LDA,N) The M-by-N matrix whose equilibration factors are to be computed. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | R | R is DOUBLE PRECISION array, dimension (M) If INFO = 0 or INFO > M, R contains the row scale factors for A. |
[out] | C | C is DOUBLE PRECISION array, dimension (N) If INFO = 0, C contains the column scale factors for A. |
[out] | ROWCND | ROWCND is DOUBLE PRECISION If INFO = 0 or INFO > M, ROWCND contains the ratio of the smallest R(i) to the largest R(i). If ROWCND >= 0.1 and AMAX is neither too large nor too small, it is not worth scaling by R. |
[out] | COLCND | COLCND is DOUBLE PRECISION If INFO = 0, COLCND contains the ratio of the smallest C(i) to the largest C(i). If COLCND >= 0.1, it is not worth scaling by C. |
[out] | AMAX | AMAX is DOUBLE PRECISION Absolute value of largest matrix element. If AMAX is very close to overflow or very close to underflow, the matrix should be scaled. |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value > 0: if INFO = i, and i is <= M: the i-th row of A is exactly zero > M: the (i-M)-th column of A is exactly zero |
Definition at line 147 of file zgeequb.f.
subroutine zgehd2 | ( | integer | N, |
integer | ILO, | ||
integer | IHI, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( * ) | TAU, | ||
complex*16, dimension( * ) | WORK, | ||
integer | INFO | ||
) |
ZGEHD2 reduces a general square matrix to upper Hessenberg form using an unblocked algorithm.
Download ZGEHD2 + dependencies [TGZ] [ZIP] [TXT]ZGEHD2 reduces a complex general matrix A to upper Hessenberg form H by a unitary similarity transformation: Q**H * A * Q = H .
[in] | N | N is INTEGER The order of the matrix A. N >= 0. |
[in] | ILO | ILO is INTEGER |
[in] | IHI | IHI is INTEGER It is assumed that A is already upper triangular in rows and columns 1:ILO-1 and IHI+1:N. ILO and IHI are normally set by a previous call to ZGEBAL; otherwise they should be set to 1 and N respectively. See Further Details. 1 <= ILO <= IHI <= max(1,N). |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the n by n general matrix to be reduced. On exit, the upper triangle and the first subdiagonal of A are overwritten with the upper Hessenberg matrix H, and the elements below the first subdiagonal, with the array TAU, represent the unitary matrix Q as a product of elementary reflectors. See Further Details. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[out] | TAU | TAU is COMPLEX*16 array, dimension (N-1) The scalar factors of the elementary reflectors (see Further Details). |
[out] | WORK | WORK is COMPLEX*16 array, dimension (N) |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value. |
The matrix Q is represented as a product of (ihi-ilo) elementary reflectors Q = H(ilo) H(ilo+1) . . . H(ihi-1). Each H(i) has the form H(i) = I - tau * v * v**H where tau is a complex scalar, and v is a complex vector with v(1:i) = 0, v(i+1) = 1 and v(ihi+1:n) = 0; v(i+2:ihi) is stored on exit in A(i+2:ihi,i), and tau in TAU(i). The contents of A are illustrated by the following example, with n = 7, ilo = 2 and ihi = 6: on entry, on exit, ( a a a a a a a ) ( a a h h h h a ) ( a a a a a a ) ( a h h h h a ) ( a a a a a a ) ( h h h h h h ) ( a a a a a a ) ( v2 h h h h h ) ( a a a a a a ) ( v2 v3 h h h h ) ( a a a a a a ) ( v2 v3 v4 h h h ) ( a ) ( a ) where a denotes an element of the original matrix A, h denotes a modified element of the upper Hessenberg matrix H, and vi denotes an element of the vector defining H(i).
Definition at line 150 of file zgehd2.f.
subroutine zgehrd | ( | integer | N, |
integer | ILO, | ||
integer | IHI, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( * ) | TAU, | ||
complex*16, dimension( * ) | WORK, | ||
integer | LWORK, | ||
integer | INFO | ||
) |
ZGEHRD
Download ZGEHRD + dependencies [TGZ] [ZIP] [TXT]ZGEHRD reduces a complex general matrix A to upper Hessenberg form H by an unitary similarity transformation: Q**H * A * Q = H .
[in] | N | N is INTEGER The order of the matrix A. N >= 0. |
[in] | ILO | ILO is INTEGER |
[in] | IHI | IHI is INTEGER It is assumed that A is already upper triangular in rows and columns 1:ILO-1 and IHI+1:N. ILO and IHI are normally set by a previous call to ZGEBAL; otherwise they should be set to 1 and N respectively. See Further Details. 1 <= ILO <= IHI <= N, if N > 0; ILO=1 and IHI=0, if N=0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the N-by-N general matrix to be reduced. On exit, the upper triangle and the first subdiagonal of A are overwritten with the upper Hessenberg matrix H, and the elements below the first subdiagonal, with the array TAU, represent the unitary matrix Q as a product of elementary reflectors. See Further Details. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[out] | TAU | TAU is COMPLEX*16 array, dimension (N-1) The scalar factors of the elementary reflectors (see Further Details). Elements 1:ILO-1 and IHI:N-1 of TAU are set to zero. |
[out] | WORK | WORK is COMPLEX*16 array, dimension (LWORK) On exit, if INFO = 0, WORK(1) returns the optimal LWORK. |
[in] | LWORK | LWORK is INTEGER The length of the array WORK. LWORK >= max(1,N). For optimum performance LWORK >= N*NB, where NB is the optimal blocksize. 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. |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value. |
The matrix Q is represented as a product of (ihi-ilo) elementary reflectors Q = H(ilo) H(ilo+1) . . . H(ihi-1). Each H(i) has the form H(i) = I - tau * v * v**H where tau is a complex scalar, and v is a complex vector with v(1:i) = 0, v(i+1) = 1 and v(ihi+1:n) = 0; v(i+2:ihi) is stored on exit in A(i+2:ihi,i), and tau in TAU(i). The contents of A are illustrated by the following example, with n = 7, ilo = 2 and ihi = 6: on entry, on exit, ( a a a a a a a ) ( a a h h h h a ) ( a a a a a a ) ( a h h h h a ) ( a a a a a a ) ( h h h h h h ) ( a a a a a a ) ( v2 h h h h h ) ( a a a a a a ) ( v2 v3 h h h h ) ( a a a a a a ) ( v2 v3 v4 h h h ) ( a ) ( a ) where a denotes an element of the original matrix A, h denotes a modified element of the upper Hessenberg matrix H, and vi denotes an element of the vector defining H(i). This file is a slight modification of LAPACK-3.0's DGEHRD subroutine incorporating improvements proposed by Quintana-Orti and Van de Geijn (2006). (See DLAHR2.)
Definition at line 169 of file zgehrd.f.
subroutine zgelq2 | ( | integer | M, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( * ) | TAU, | ||
complex*16, dimension( * ) | WORK, | ||
integer | INFO | ||
) |
ZGELQ2 computes the LQ factorization of a general rectangular matrix using an unblocked algorithm.
Download ZGELQ2 + dependencies [TGZ] [ZIP] [TXT]ZGELQ2 computes an LQ factorization of a complex m by n matrix A: A = L * Q.
[in] | M | M is INTEGER The number of rows of the matrix A. M >= 0. |
[in] | N | N is INTEGER The number of columns of the matrix A. N >= 0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the m by n matrix A. On exit, the elements on and below the diagonal of the array contain the m by min(m,n) lower trapezoidal matrix L (L is lower triangular if m <= n); the elements above the diagonal, with the array TAU, represent the unitary matrix Q as a product of elementary reflectors (see Further Details). |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | TAU | TAU is COMPLEX*16 array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details). |
[out] | WORK | WORK is COMPLEX*16 array, dimension (M) |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value |
The matrix Q is represented as a product of elementary reflectors Q = H(k)**H . . . H(2)**H H(1)**H, where k = min(m,n). Each H(i) has the form H(i) = I - tau * v * v**H where tau is a complex scalar, and v is a complex vector with v(1:i-1) = 0 and v(i) = 1; conjg(v(i+1:n)) is stored on exit in A(i,i+1:n), and tau in TAU(i).
Definition at line 122 of file zgelq2.f.
subroutine zgelqf | ( | integer | M, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( * ) | TAU, | ||
complex*16, dimension( * ) | WORK, | ||
integer | LWORK, | ||
integer | INFO | ||
) |
ZGELQF
Download ZGELQF + dependencies [TGZ] [ZIP] [TXT]ZGELQF computes an LQ factorization of a complex M-by-N matrix A: A = L * Q.
[in] | M | M is INTEGER The number of rows of the matrix A. M >= 0. |
[in] | N | N is INTEGER The number of columns of the matrix A. N >= 0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the M-by-N matrix A. On exit, the elements on and below the diagonal of the array contain the m-by-min(m,n) lower trapezoidal matrix L (L is lower triangular if m <= n); the elements above the diagonal, with the array TAU, represent the unitary matrix Q as a product of elementary reflectors (see Further Details). |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | TAU | TAU is COMPLEX*16 array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details). |
[out] | WORK | WORK is COMPLEX*16 array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK(1) returns the optimal LWORK. |
[in] | LWORK | LWORK is INTEGER The dimension of the array WORK. LWORK >= max(1,M). For optimum performance LWORK >= M*NB, where NB is the optimal blocksize. 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. |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value |
The matrix Q is represented as a product of elementary reflectors Q = H(k)**H . . . H(2)**H H(1)**H, where k = min(m,n). Each H(i) has the form H(i) = I - tau * v * v**H where tau is a complex scalar, and v is a complex vector with v(1:i-1) = 0 and v(i) = 1; conjg(v(i+1:n)) is stored on exit in A(i,i+1:n), and tau in TAU(i).
Definition at line 136 of file zgelqf.f.
subroutine zgemqrt | ( | character | SIDE, |
character | TRANS, | ||
integer | M, | ||
integer | N, | ||
integer | K, | ||
integer | NB, | ||
complex*16, dimension( ldv, * ) | V, | ||
integer | LDV, | ||
complex*16, dimension( ldt, * ) | T, | ||
integer | LDT, | ||
complex*16, dimension( ldc, * ) | C, | ||
integer | LDC, | ||
complex*16, dimension( * ) | WORK, | ||
integer | INFO | ||
) |
ZGEMQRT
Download ZGEMQRT + dependencies [TGZ] [ZIP] [TXT]ZGEMQRT overwrites the general complex M-by-N matrix C with SIDE = 'L' SIDE = 'R' TRANS = 'N': Q C C Q TRANS = 'C': Q**H C C Q**H where Q is a complex orthogonal matrix defined as the product of K elementary reflectors: Q = H(1) H(2) . . . H(K) = I - V T V**H generated using the compact WY representation as returned by ZGEQRT. Q is of order M if SIDE = 'L' and of order N if SIDE = 'R'.
[in] | SIDE | SIDE is CHARACTER*1 = 'L': apply Q or Q**H from the Left; = 'R': apply Q or Q**H from the Right. |
[in] | TRANS | TRANS is CHARACTER*1 = 'N': No transpose, apply Q; = 'C': Transpose, apply Q**H. |
[in] | M | M is INTEGER The number of rows of the matrix C. M >= 0. |
[in] | N | N is INTEGER The number of columns of the matrix C. N >= 0. |
[in] | K | K is INTEGER The number of elementary reflectors whose product defines the matrix Q. If SIDE = 'L', M >= K >= 0; if SIDE = 'R', N >= K >= 0. |
[in] | NB | NB is INTEGER The block size used for the storage of T. K >= NB >= 1. This must be the same value of NB used to generate T in CGEQRT. |
[in] | V | V is COMPLEX*16 array, dimension (LDV,K) The i-th column must contain the vector which defines the elementary reflector H(i), for i = 1,2,...,k, as returned by CGEQRT in the first K columns of its array argument A. |
[in] | LDV | LDV is INTEGER The leading dimension of the array V. If SIDE = 'L', LDA >= max(1,M); if SIDE = 'R', LDA >= max(1,N). |
[in] | T | T is COMPLEX*16 array, dimension (LDT,K) The upper triangular factors of the block reflectors as returned by CGEQRT, stored as a NB-by-N matrix. |
[in] | LDT | LDT is INTEGER The leading dimension of the array T. LDT >= NB. |
[in,out] | C | C is COMPLEX*16 array, dimension (LDC,N) On entry, the M-by-N matrix C. On exit, C is overwritten by Q C, Q**H C, C Q**H or C Q. |
[in] | LDC | LDC is INTEGER The leading dimension of the array C. LDC >= max(1,M). |
[out] | WORK | WORK is COMPLEX*16 array. The dimension of WORK is N*NB if SIDE = 'L', or M*NB if SIDE = 'R'. |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value |
Definition at line 168 of file zgemqrt.f.
subroutine zgeql2 | ( | integer | M, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( * ) | TAU, | ||
complex*16, dimension( * ) | WORK, | ||
integer | INFO | ||
) |
ZGEQL2 computes the QL factorization of a general rectangular matrix using an unblocked algorithm.
Download ZGEQL2 + dependencies [TGZ] [ZIP] [TXT]ZGEQL2 computes a QL factorization of a complex m by n matrix A: A = Q * L.
[in] | M | M is INTEGER The number of rows of the matrix A. M >= 0. |
[in] | N | N is INTEGER The number of columns of the matrix A. N >= 0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the m by n matrix A. On exit, if m >= n, the lower triangle of the subarray A(m-n+1:m,1:n) contains the n by n lower triangular matrix L; if m <= n, the elements on and below the (n-m)-th superdiagonal contain the m by n lower trapezoidal matrix L; the remaining elements, with the array TAU, represent the unitary matrix Q as a product of elementary reflectors (see Further Details). |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | TAU | TAU is COMPLEX*16 array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details). |
[out] | WORK | WORK is COMPLEX*16 array, dimension (N) |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value |
The matrix Q is represented as a product of elementary reflectors Q = H(k) . . . H(2) H(1), where k = min(m,n). Each H(i) has the form H(i) = I - tau * v * v**H where tau is a complex scalar, and v is a complex vector with v(m-k+i+1:m) = 0 and v(m-k+i) = 1; v(1:m-k+i-1) is stored on exit in A(1:m-k+i-1,n-k+i), and tau in TAU(i).
Definition at line 124 of file zgeql2.f.
subroutine zgeqlf | ( | integer | M, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( * ) | TAU, | ||
complex*16, dimension( * ) | WORK, | ||
integer | LWORK, | ||
integer | INFO | ||
) |
ZGEQLF
Download ZGEQLF + dependencies [TGZ] [ZIP] [TXT]ZGEQLF computes a QL factorization of a complex M-by-N matrix A: A = Q * L.
[in] | M | M is INTEGER The number of rows of the matrix A. M >= 0. |
[in] | N | N is INTEGER The number of columns of the matrix A. N >= 0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the M-by-N matrix A. On exit, if m >= n, the lower triangle of the subarray A(m-n+1:m,1:n) contains the N-by-N lower triangular matrix L; if m <= n, the elements on and below the (n-m)-th superdiagonal contain the M-by-N lower trapezoidal matrix L; the remaining elements, with the array TAU, represent the unitary matrix Q as a product of elementary reflectors (see Further Details). |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | TAU | TAU is COMPLEX*16 array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details). |
[out] | WORK | WORK is COMPLEX*16 array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK(1) returns the optimal LWORK. |
[in] | LWORK | LWORK is INTEGER The dimension of the array WORK. LWORK >= max(1,N). For optimum performance LWORK >= N*NB, where NB is the optimal blocksize. 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. |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value |
The matrix Q is represented as a product of elementary reflectors Q = H(k) . . . H(2) H(1), where k = min(m,n). Each H(i) has the form H(i) = I - tau * v * v**H where tau is a complex scalar, and v is a complex vector with v(m-k+i+1:m) = 0 and v(m-k+i) = 1; v(1:m-k+i-1) is stored on exit in A(1:m-k+i-1,n-k+i), and tau in TAU(i).
Definition at line 139 of file zgeqlf.f.
subroutine zgeqp3 | ( | integer | M, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
integer, dimension( * ) | JPVT, | ||
complex*16, dimension( * ) | TAU, | ||
complex*16, dimension( * ) | WORK, | ||
integer | LWORK, | ||
double precision, dimension( * ) | RWORK, | ||
integer | INFO | ||
) |
ZGEQP3
Download ZGEQP3 + dependencies [TGZ] [ZIP] [TXT]ZGEQP3 computes a QR factorization with column pivoting of a matrix A: A*P = Q*R using Level 3 BLAS.
[in] | M | M is INTEGER The number of rows of the matrix A. M >= 0. |
[in] | N | N is INTEGER The number of columns of the matrix A. N >= 0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the M-by-N matrix A. On exit, the upper triangle of the array contains the min(M,N)-by-N upper trapezoidal matrix R; the elements below the diagonal, together with the array TAU, represent the unitary matrix Q as a product of min(M,N) elementary reflectors. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[in,out] | JPVT | JPVT is INTEGER array, dimension (N) On entry, if JPVT(J).ne.0, the J-th column of A is permuted to the front of A*P (a leading column); if JPVT(J)=0, the J-th column of A is a free column. On exit, if JPVT(J)=K, then the J-th column of A*P was the the K-th column of A. |
[out] | TAU | TAU is COMPLEX*16 array, dimension (min(M,N)) The scalar factors of the elementary reflectors. |
[out] | WORK | WORK is COMPLEX*16 array, dimension (MAX(1,LWORK)) On exit, if INFO=0, WORK(1) returns the optimal LWORK. |
[in] | LWORK | LWORK is INTEGER The dimension of the array WORK. LWORK >= N+1. For optimal performance LWORK >= ( N+1 )*NB, where NB is the optimal blocksize. 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. |
[out] | RWORK | RWORK is DOUBLE PRECISION array, dimension (2*N) |
[out] | INFO | INFO is INTEGER = 0: successful exit. < 0: if INFO = -i, the i-th argument had an illegal value. |
The matrix Q is represented as a product of elementary reflectors Q = H(1) H(2) . . . H(k), where k = min(m,n). Each H(i) has the form H(i) = I - tau * v * v**H where tau is a complex scalar, and v is a real/complex vector with v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i), and tau in TAU(i).
Definition at line 159 of file zgeqp3.f.
subroutine zgeqpf | ( | integer | M, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
integer, dimension( * ) | JPVT, | ||
complex*16, dimension( * ) | TAU, | ||
complex*16, dimension( * ) | WORK, | ||
double precision, dimension( * ) | RWORK, | ||
integer | INFO | ||
) |
ZGEQPF
Download ZGEQPF + dependencies [TGZ] [ZIP] [TXT]This routine is deprecated and has been replaced by routine ZGEQP3. ZGEQPF computes a QR factorization with column pivoting of a complex M-by-N matrix A: A*P = Q*R.
[in] | M | M is INTEGER The number of rows of the matrix A. M >= 0. |
[in] | N | N is INTEGER The number of columns of the matrix A. N >= 0 |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the M-by-N matrix A. On exit, the upper triangle of the array contains the min(M,N)-by-N upper triangular matrix R; the elements below the diagonal, together with the array TAU, represent the unitary matrix Q as a product of min(m,n) elementary reflectors. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[in,out] | JPVT | JPVT is INTEGER array, dimension (N) On entry, if JPVT(i) .ne. 0, the i-th column of A is permuted to the front of A*P (a leading column); if JPVT(i) = 0, the i-th column of A is a free column. On exit, if JPVT(i) = k, then the i-th column of A*P was the k-th column of A. |
[out] | TAU | TAU is COMPLEX*16 array, dimension (min(M,N)) The scalar factors of the elementary reflectors. |
[out] | WORK | WORK is COMPLEX*16 array, dimension (N) |
[out] | RWORK | RWORK is DOUBLE PRECISION array, dimension (2*N) |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value |
The matrix Q is represented as a product of elementary reflectors Q = H(1) H(2) . . . H(n) Each H(i) has the form H = I - tau * v * v**H where tau is a complex scalar, and v is a complex vector with v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i). The matrix P is represented in jpvt as follows: If jpvt(j) = i then the jth column of P is the ith canonical unit vector. Partial column norm updating strategy modified by Z. Drmac and Z. Bujanovic, Dept. of Mathematics, University of Zagreb, Croatia. -- April 2011 -- For more details see LAPACK Working Note 176.
Definition at line 149 of file zgeqpf.f.
subroutine zgeqr2 | ( | integer | M, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( * ) | TAU, | ||
complex*16, dimension( * ) | WORK, | ||
integer | INFO | ||
) |
ZGEQR2 computes the QR factorization of a general rectangular matrix using an unblocked algorithm.
Download ZGEQR2 + dependencies [TGZ] [ZIP] [TXT]ZGEQR2 computes a QR factorization of a complex m by n matrix A: A = Q * R.
[in] | M | M is INTEGER The number of rows of the matrix A. M >= 0. |
[in] | N | N is INTEGER The number of columns of the matrix A. N >= 0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the m by n matrix A. On exit, the elements on and above the diagonal of the array contain the min(m,n) by n upper trapezoidal matrix R (R is upper triangular if m >= n); the elements below the diagonal, with the array TAU, represent the unitary matrix Q as a product of elementary reflectors (see Further Details). |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | TAU | TAU is COMPLEX*16 array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details). |
[out] | WORK | WORK is COMPLEX*16 array, dimension (N) |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value |
The matrix Q is represented as a product of elementary reflectors Q = H(1) H(2) . . . H(k), where k = min(m,n). Each H(i) has the form H(i) = I - tau * v * v**H where tau is a complex scalar, and v is a complex vector with v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i), and tau in TAU(i).
Definition at line 122 of file zgeqr2.f.
subroutine zgeqr2p | ( | integer | M, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( * ) | TAU, | ||
complex*16, dimension( * ) | WORK, | ||
integer | INFO | ||
) |
ZGEQR2P computes the QR factorization of a general rectangular matrix with non-negative diagonal elements using an unblocked algorithm.
Download ZGEQR2P + dependencies [TGZ] [ZIP] [TXT]ZGEQR2P computes a QR factorization of a complex m by n matrix A: A = Q * R.
[in] | M | M is INTEGER The number of rows of the matrix A. M >= 0. |
[in] | N | N is INTEGER The number of columns of the matrix A. N >= 0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the m by n matrix A. On exit, the elements on and above the diagonal of the array contain the min(m,n) by n upper trapezoidal matrix R (R is upper triangular if m >= n); the elements below the diagonal, with the array TAU, represent the unitary matrix Q as a product of elementary reflectors (see Further Details). |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | TAU | TAU is COMPLEX*16 array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details). |
[out] | WORK | WORK is COMPLEX*16 array, dimension (N) |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value |
The matrix Q is represented as a product of elementary reflectors Q = H(1) H(2) . . . H(k), where k = min(m,n). Each H(i) has the form H(i) = I - tau * v * v**H where tau is a complex scalar, and v is a complex vector with v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i), and tau in TAU(i).
Definition at line 122 of file zgeqr2p.f.
subroutine zgeqrf | ( | integer | M, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( * ) | TAU, | ||
complex*16, dimension( * ) | WORK, | ||
integer | LWORK, | ||
integer | INFO | ||
) |
ZGEQRF
Download ZGEQRF + dependencies [TGZ] [ZIP] [TXT]ZGEQRF computes a QR factorization of a complex M-by-N matrix A: A = Q * R.
[in] | M | M is INTEGER The number of rows of the matrix A. M >= 0. |
[in] | N | N is INTEGER The number of columns of the matrix A. N >= 0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the M-by-N matrix A. On exit, the elements on and above the diagonal of the array contain the min(M,N)-by-N upper trapezoidal matrix R (R is upper triangular if m >= n); the elements below the diagonal, with the array TAU, represent the unitary matrix Q as a product of min(m,n) elementary reflectors (see Further Details). |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | TAU | TAU is COMPLEX*16 array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details). |
[out] | WORK | WORK is COMPLEX*16 array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK(1) returns the optimal LWORK. |
[in] | LWORK | LWORK is INTEGER The dimension of the array WORK. LWORK >= max(1,N). For optimum performance LWORK >= N*NB, where NB is the optimal blocksize. 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. |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value |
The matrix Q is represented as a product of elementary reflectors Q = H(1) H(2) . . . H(k), where k = min(m,n). Each H(i) has the form H(i) = I - tau * v * v**H where tau is a complex scalar, and v is a complex vector with v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i), and tau in TAU(i).
Definition at line 137 of file zgeqrf.f.
subroutine zgeqrfp | ( | integer | M, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( * ) | TAU, | ||
complex*16, dimension( * ) | WORK, | ||
integer | LWORK, | ||
integer | INFO | ||
) |
ZGEQRFP
Download ZGEQRFP + dependencies [TGZ] [ZIP] [TXT]ZGEQRFP computes a QR factorization of a complex M-by-N matrix A: A = Q * R.
[in] | M | M is INTEGER The number of rows of the matrix A. M >= 0. |
[in] | N | N is INTEGER The number of columns of the matrix A. N >= 0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the M-by-N matrix A. On exit, the elements on and above the diagonal of the array contain the min(M,N)-by-N upper trapezoidal matrix R (R is upper triangular if m >= n); the elements below the diagonal, with the array TAU, represent the unitary matrix Q as a product of min(m,n) elementary reflectors (see Further Details). |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | TAU | TAU is COMPLEX*16 array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details). |
[out] | WORK | WORK is COMPLEX*16 array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK(1) returns the optimal LWORK. |
[in] | LWORK | LWORK is INTEGER The dimension of the array WORK. LWORK >= max(1,N). For optimum performance LWORK >= N*NB, where NB is the optimal blocksize. 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. |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value |
The matrix Q is represented as a product of elementary reflectors Q = H(1) H(2) . . . H(k), where k = min(m,n). Each H(i) has the form H(i) = I - tau * v * v**H where tau is a complex scalar, and v is a complex vector with v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i), and tau in TAU(i).
Definition at line 137 of file zgeqrfp.f.
subroutine zgeqrt | ( | integer | M, |
integer | N, | ||
integer | NB, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( ldt, * ) | T, | ||
integer | LDT, | ||
complex*16, dimension( * ) | WORK, | ||
integer | INFO | ||
) |
ZGEQRT
Download ZGEQRT + dependencies [TGZ] [ZIP] [TXT]ZGEQRT computes a blocked QR factorization of a complex M-by-N matrix A using the compact WY representation of Q.
[in] | M | M is INTEGER The number of rows of the matrix A. M >= 0. |
[in] | N | N is INTEGER The number of columns of the matrix A. N >= 0. |
[in] | NB | NB is INTEGER The block size to be used in the blocked QR. MIN(M,N) >= NB >= 1. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the M-by-N matrix A. On exit, the elements on and above the diagonal of the array contain the min(M,N)-by-N upper trapezoidal matrix R (R is upper triangular if M >= N); the elements below the diagonal are the columns of V. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | T | T is COMPLEX*16 array, dimension (LDT,MIN(M,N)) The upper triangular block reflectors stored in compact form as a sequence of upper triangular blocks. See below for further details. |
[in] | LDT | LDT is INTEGER The leading dimension of the array T. LDT >= NB. |
[out] | WORK | WORK is COMPLEX*16 array, dimension (NB*N) |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value |
The matrix V stores the elementary reflectors H(i) in the i-th column below the diagonal. For example, if M=5 and N=3, the matrix V is V = ( 1 ) ( v1 1 ) ( v1 v2 1 ) ( v1 v2 v3 ) ( v1 v2 v3 ) where the vi's represent the vectors which define H(i), which are returned in the matrix A. The 1's along the diagonal of V are not stored in A. Let K=MIN(M,N). The number of blocks is B = ceiling(K/NB), where each block is of order NB except for the last block, which is of order IB = K - (B-1)*NB. For each of the B blocks, a upper triangular block reflector factor is computed: T1, T2, ..., TB. The NB-by-NB (and IB-by-IB for the last block) T's are stored in the NB-by-N matrix T as T = (T1 T2 ... TB).
Definition at line 142 of file zgeqrt.f.
subroutine zgeqrt2 | ( | integer | M, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( ldt, * ) | T, | ||
integer | LDT, | ||
integer | INFO | ||
) |
ZGEQRT2 computes a QR factorization of a general real or complex matrix using the compact WY representation of Q.
Download ZGEQRT2 + dependencies [TGZ] [ZIP] [TXT]ZGEQRT2 computes a QR factorization of a complex M-by-N matrix A, using the compact WY representation of Q.
[in] | M | M is INTEGER The number of rows of the matrix A. M >= N. |
[in] | N | N is INTEGER The number of columns of the matrix A. N >= 0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the complex M-by-N matrix A. On exit, the elements on and above the diagonal contain the N-by-N upper triangular matrix R; the elements below the diagonal are the columns of V. See below for further details. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | T | T is COMPLEX*16 array, dimension (LDT,N) The N-by-N upper triangular factor of the block reflector. The elements on and above the diagonal contain the block reflector T; the elements below the diagonal are not used. See below for further details. |
[in] | LDT | LDT is INTEGER The leading dimension of the array T. LDT >= max(1,N). |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value |
The matrix V stores the elementary reflectors H(i) in the i-th column below the diagonal. For example, if M=5 and N=3, the matrix V is V = ( 1 ) ( v1 1 ) ( v1 v2 1 ) ( v1 v2 v3 ) ( v1 v2 v3 ) where the vi's represent the vectors which define H(i), which are returned in the matrix A. The 1's along the diagonal of V are not stored in A. The block reflector H is then given by H = I - V * T * V**H where V**H is the conjugate transpose of V.
Definition at line 128 of file zgeqrt2.f.
recursive subroutine zgeqrt3 | ( | integer | M, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( ldt, * ) | T, | ||
integer | LDT, | ||
integer | INFO | ||
) |
ZGEQRT3 recursively computes a QR factorization of a general real or complex matrix using the compact WY representation of Q.
Download ZGEQRT3 + dependencies [TGZ] [ZIP] [TXT]ZGEQRT3 recursively computes a QR factorization of a complex M-by-N matrix A, using the compact WY representation of Q. Based on the algorithm of Elmroth and Gustavson, IBM J. Res. Develop. Vol 44 No. 4 July 2000.
[in] | M | M is INTEGER The number of rows of the matrix A. M >= N. |
[in] | N | N is INTEGER The number of columns of the matrix A. N >= 0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the complex M-by-N matrix A. On exit, the elements on and above the diagonal contain the N-by-N upper triangular matrix R; the elements below the diagonal are the columns of V. See below for further details. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | T | T is COMPLEX*16 array, dimension (LDT,N) The N-by-N upper triangular factor of the block reflector. The elements on and above the diagonal contain the block reflector T; the elements below the diagonal are not used. See below for further details. |
[in] | LDT | LDT is INTEGER The leading dimension of the array T. LDT >= max(1,N). |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value |
The matrix V stores the elementary reflectors H(i) in the i-th column below the diagonal. For example, if M=5 and N=3, the matrix V is V = ( 1 ) ( v1 1 ) ( v1 v2 1 ) ( v1 v2 v3 ) ( v1 v2 v3 ) where the vi's represent the vectors which define H(i), which are returned in the matrix A. The 1's along the diagonal of V are not stored in A. The block reflector H is then given by H = I - V * T * V**H where V**H is the conjugate transpose of V. For details of the algorithm, see Elmroth and Gustavson (cited above).
Definition at line 133 of file zgeqrt3.f.
subroutine zgerfs | ( | character | TRANS, |
integer | N, | ||
integer | NRHS, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( ldaf, * ) | AF, | ||
integer | LDAF, | ||
integer, dimension( * ) | IPIV, | ||
complex*16, dimension( ldb, * ) | B, | ||
integer | LDB, | ||
complex*16, dimension( ldx, * ) | X, | ||
integer | LDX, | ||
double precision, dimension( * ) | FERR, | ||
double precision, dimension( * ) | BERR, | ||
complex*16, dimension( * ) | WORK, | ||
double precision, dimension( * ) | RWORK, | ||
integer | INFO | ||
) |
ZGERFS
Download ZGERFS + dependencies [TGZ] [ZIP] [TXT]ZGERFS improves the computed solution to a system of linear equations and provides error bounds and backward error estimates for the solution.
[in] | TRANS | TRANS is 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) |
[in] | N | N is INTEGER The order of the matrix A. N >= 0. |
[in] | NRHS | NRHS is INTEGER The number of right hand sides, i.e., the number of columns of the matrices B and X. NRHS >= 0. |
[in] | A | A is COMPLEX*16 array, dimension (LDA,N) The original N-by-N matrix A. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in] | AF | AF is COMPLEX*16 array, dimension (LDAF,N) The factors L and U from the factorization A = P*L*U as computed by ZGETRF. |
[in] | LDAF | LDAF is INTEGER The leading dimension of the array AF. LDAF >= max(1,N). |
[in] | IPIV | IPIV is INTEGER array, dimension (N) The pivot indices from ZGETRF; for 1<=i<=N, row i of the matrix was interchanged with row IPIV(i). |
[in] | B | B is COMPLEX*16 array, dimension (LDB,NRHS) The right hand side matrix B. |
[in] | LDB | LDB is INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in,out] | X | X is COMPLEX*16 array, dimension (LDX,NRHS) On entry, the solution matrix X, as computed by ZGETRS. On exit, the improved solution matrix X. |
[in] | LDX | LDX is INTEGER The leading dimension of the array X. LDX >= max(1,N). |
[out] | FERR | FERR is 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. |
[out] | BERR | BERR is 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). |
[out] | WORK | WORK is COMPLEX*16 array, dimension (2*N) |
[out] | RWORK | RWORK is DOUBLE PRECISION array, dimension (N) |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value |
ITMAX is the maximum number of steps of iterative refinement.
Definition at line 186 of file zgerfs.f.
subroutine zgerfsx | ( | character | TRANS, |
character | EQUED, | ||
integer | N, | ||
integer | NRHS, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( ldaf, * ) | AF, | ||
integer | LDAF, | ||
integer, dimension( * ) | IPIV, | ||
double precision, dimension( * ) | R, | ||
double precision, dimension( * ) | C, | ||
complex*16, dimension( ldb, * ) | B, | ||
integer | LDB, | ||
complex*16, dimension( ldx , * ) | X, | ||
integer | LDX, | ||
double precision | RCOND, | ||
double precision, dimension( * ) | BERR, | ||
integer | N_ERR_BNDS, | ||
double precision, dimension( nrhs, * ) | ERR_BNDS_NORM, | ||
double precision, dimension( nrhs, * ) | ERR_BNDS_COMP, | ||
integer | NPARAMS, | ||
double precision, dimension( * ) | PARAMS, | ||
complex*16, dimension( * ) | WORK, | ||
double precision, dimension( * ) | RWORK, | ||
integer | INFO | ||
) |
ZGERFSX
Download ZGERFSX + dependencies [TGZ] [ZIP] [TXT]ZGERFSX improves the computed solution to a system of linear equations and provides error bounds and backward error estimates for the solution. In addition to normwise error bound, the code provides maximum componentwise error bound if possible. See comments for ERR_BNDS_NORM and ERR_BNDS_COMP for details of the error bounds. The original system of linear equations may have been equilibrated before calling this routine, as described by arguments EQUED, R and C below. In this case, the solution and error bounds returned are for the original unequilibrated system.
Some optional parameters are bundled in the PARAMS array. These settings determine how refinement is performed, but often the defaults are acceptable. If the defaults are acceptable, users can pass NPARAMS = 0 which prevents the source code from accessing the PARAMS argument.
[in] | TRANS | TRANS is 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) |
[in] | EQUED | EQUED is CHARACTER*1 Specifies the form of equilibration that was done to A before calling this routine. This is needed to compute the solution and error bounds correctly. = 'N': No equilibration = 'R': Row equilibration, i.e., A has been premultiplied by diag(R). = 'C': Column equilibration, i.e., A has been postmultiplied by diag(C). = 'B': Both row and column equilibration, i.e., A has been replaced by diag(R) * A * diag(C). The right hand side B has been changed accordingly. |
[in] | N | N is INTEGER The order of the matrix A. N >= 0. |
[in] | NRHS | NRHS is INTEGER The number of right hand sides, i.e., the number of columns of the matrices B and X. NRHS >= 0. |
[in] | A | A is COMPLEX*16 array, dimension (LDA,N) The original N-by-N matrix A. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in] | AF | AF is COMPLEX*16 array, dimension (LDAF,N) The factors L and U from the factorization A = P*L*U as computed by ZGETRF. |
[in] | LDAF | LDAF is INTEGER The leading dimension of the array AF. LDAF >= max(1,N). |
[in] | IPIV | IPIV is INTEGER array, dimension (N) The pivot indices from ZGETRF; for 1<=i<=N, row i of the matrix was interchanged with row IPIV(i). |
[in] | R | R is DOUBLE PRECISION array, dimension (N) The row scale factors for A. If EQUED = 'R' or 'B', A is multiplied on the left by diag(R); if EQUED = 'N' or 'C', R is not accessed. If R is accessed, each element of R should be a power of the radix to ensure a reliable solution and error estimates. Scaling by powers of the radix does not cause rounding errors unless the result underflows or overflows. Rounding errors during scaling lead to refining with a matrix that is not equivalent to the input matrix, producing error estimates that may not be reliable. |
[in] | C | C is DOUBLE PRECISION array, dimension (N) The column scale factors for A. If EQUED = 'C' or 'B', A is multiplied on the right by diag(C); if EQUED = 'N' or 'R', C is not accessed. If C is accessed, each element of C should be a power of the radix to ensure a reliable solution and error estimates. Scaling by powers of the radix does not cause rounding errors unless the result underflows or overflows. Rounding errors during scaling lead to refining with a matrix that is not equivalent to the input matrix, producing error estimates that may not be reliable. |
[in] | B | B is COMPLEX*16 array, dimension (LDB,NRHS) The right hand side matrix B. |
[in] | LDB | LDB is INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in,out] | X | X is COMPLEX*16 array, dimension (LDX,NRHS) On entry, the solution matrix X, as computed by ZGETRS. On exit, the improved solution matrix X. |
[in] | LDX | LDX is INTEGER The leading dimension of the array X. LDX >= max(1,N). |
[out] | RCOND | RCOND is DOUBLE PRECISION Reciprocal scaled condition number. This is an estimate of the reciprocal Skeel condition number of the matrix A after equilibration (if done). If this is less than the machine precision (in particular, if it is zero), the matrix is singular to working precision. Note that the error may still be small even if this number is very small and the matrix appears ill- conditioned. |
[out] | BERR | BERR is DOUBLE PRECISION array, dimension (NRHS) Componentwise relative backward error. This is 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). |
[in] | N_ERR_BNDS | N_ERR_BNDS is INTEGER Number of error bounds to return for each right hand side and each type (normwise or componentwise). See ERR_BNDS_NORM and ERR_BNDS_COMP below. |
[out] | ERR_BNDS_NORM | ERR_BNDS_NORM is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS) For each right-hand side, this array contains information about various error bounds and condition numbers corresponding to the normwise relative error, which is defined as follows: Normwise relative error in the ith solution vector: max_j (abs(XTRUE(j,i) - X(j,i))) ------------------------------ max_j abs(X(j,i)) The array is indexed by the type of error information as described below. There currently are up to three pieces of information returned. The first index in ERR_BNDS_NORM(i,:) corresponds to the ith right-hand side. The second index in ERR_BNDS_NORM(:,err) contains the following three fields: err = 1 "Trust/don't trust" boolean. Trust the answer if the reciprocal condition number is less than the threshold sqrt(n) * dlamch('Epsilon'). err = 2 "Guaranteed" error bound: The estimated forward error, almost certainly within a factor of 10 of the true error so long as the next entry is greater than the threshold sqrt(n) * dlamch('Epsilon'). This error bound should only be trusted if the previous boolean is true. err = 3 Reciprocal condition number: Estimated normwise reciprocal condition number. Compared with the threshold sqrt(n) * dlamch('Epsilon') to determine if the error estimate is "guaranteed". These reciprocal condition numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some appropriately scaled matrix Z. Let Z = S*A, where S scales each row by a power of the radix so all absolute row sums of Z are approximately 1. See Lapack Working Note 165 for further details and extra cautions. |
[out] | ERR_BNDS_COMP | ERR_BNDS_COMP is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS) For each right-hand side, this array contains information about various error bounds and condition numbers corresponding to the componentwise relative error, which is defined as follows: Componentwise relative error in the ith solution vector: abs(XTRUE(j,i) - X(j,i)) max_j ---------------------- abs(X(j,i)) The array is indexed by the right-hand side i (on which the componentwise relative error depends), and the type of error information as described below. There currently are up to three pieces of information returned for each right-hand side. If componentwise accuracy is not requested (PARAMS(3) = 0.0), then ERR_BNDS_COMP is not accessed. If N_ERR_BNDS .LT. 3, then at most the first (:,N_ERR_BNDS) entries are returned. The first index in ERR_BNDS_COMP(i,:) corresponds to the ith right-hand side. The second index in ERR_BNDS_COMP(:,err) contains the following three fields: err = 1 "Trust/don't trust" boolean. Trust the answer if the reciprocal condition number is less than the threshold sqrt(n) * dlamch('Epsilon'). err = 2 "Guaranteed" error bound: The estimated forward error, almost certainly within a factor of 10 of the true error so long as the next entry is greater than the threshold sqrt(n) * dlamch('Epsilon'). This error bound should only be trusted if the previous boolean is true. err = 3 Reciprocal condition number: Estimated componentwise reciprocal condition number. Compared with the threshold sqrt(n) * dlamch('Epsilon') to determine if the error estimate is "guaranteed". These reciprocal condition numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some appropriately scaled matrix Z. Let Z = S*(A*diag(x)), where x is the solution for the current right-hand side and S scales each row of A*diag(x) by a power of the radix so all absolute row sums of Z are approximately 1. See Lapack Working Note 165 for further details and extra cautions. |
[in] | NPARAMS | NPARAMS is INTEGER Specifies the number of parameters set in PARAMS. If .LE. 0, the PARAMS array is never referenced and default values are used. |
[in,out] | PARAMS | PARAMS is DOUBLE PRECISION array, dimension NPARAMS Specifies algorithm parameters. If an entry is .LT. 0.0, then that entry will be filled with default value used for that parameter. Only positions up to NPARAMS are accessed; defaults are used for higher-numbered parameters. PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative refinement or not. Default: 1.0D+0 = 0.0 : No refinement is performed, and no error bounds are computed. = 1.0 : Use the double-precision refinement algorithm, possibly with doubled-single computations if the compilation environment does not support DOUBLE PRECISION. (other values are reserved for future use) PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual computations allowed for refinement. Default: 10 Aggressive: Set to 100 to permit convergence using approximate factorizations or factorizations other than LU. If the factorization uses a technique other than Gaussian elimination, the guarantees in err_bnds_norm and err_bnds_comp may no longer be trustworthy. PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code will attempt to find a solution with small componentwise relative error in the double-precision algorithm. Positive is true, 0.0 is false. Default: 1.0 (attempt componentwise convergence) |
[out] | WORK | WORK is COMPLEX*16 array, dimension (2*N) |
[out] | RWORK | RWORK is DOUBLE PRECISION array, dimension (2*N) |
[out] | INFO | INFO is INTEGER = 0: Successful exit. The solution to every right-hand side is guaranteed. < 0: If INFO = -i, the i-th argument had an illegal value > 0 and <= N: U(INFO,INFO) is exactly zero. The factorization has been completed, but the factor U is exactly singular, so the solution and error bounds could not be computed. RCOND = 0 is returned. = N+J: The solution corresponding to the Jth right-hand side is not guaranteed. The solutions corresponding to other right- hand sides K with K > J may not be guaranteed as well, but only the first such right-hand side is reported. If a small componentwise error is not requested (PARAMS(3) = 0.0) then the Jth right-hand side is the first with a normwise error bound that is not guaranteed (the smallest J such that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0) the Jth right-hand side is the first with either a normwise or componentwise error bound that is not guaranteed (the smallest J such that either ERR_BNDS_NORM(J,1) = 0.0 or ERR_BNDS_COMP(J,1) = 0.0). See the definition of ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information about all of the right-hand sides check ERR_BNDS_NORM or ERR_BNDS_COMP. |
Definition at line 412 of file zgerfsx.f.
subroutine zgerq2 | ( | integer | M, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( * ) | TAU, | ||
complex*16, dimension( * ) | WORK, | ||
integer | INFO | ||
) |
ZGERQ2 computes the RQ factorization of a general rectangular matrix using an unblocked algorithm.
Download ZGERQ2 + dependencies [TGZ] [ZIP] [TXT]ZGERQ2 computes an RQ factorization of a complex m by n matrix A: A = R * Q.
[in] | M | M is INTEGER The number of rows of the matrix A. M >= 0. |
[in] | N | N is INTEGER The number of columns of the matrix A. N >= 0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the m by n matrix A. On exit, if m <= n, the upper triangle of the subarray A(1:m,n-m+1:n) contains the m by m upper triangular matrix R; if m >= n, the elements on and above the (m-n)-th subdiagonal contain the m by n upper trapezoidal matrix R; the remaining elements, with the array TAU, represent the unitary matrix Q as a product of elementary reflectors (see Further Details). |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | TAU | TAU is COMPLEX*16 array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details). |
[out] | WORK | WORK is COMPLEX*16 array, dimension (M) |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value |
The matrix Q is represented as a product of elementary reflectors Q = H(1)**H H(2)**H . . . H(k)**H, where k = min(m,n). Each H(i) has the form H(i) = I - tau * v * v**H where tau is a complex scalar, and v is a complex vector with v(n-k+i+1:n) = 0 and v(n-k+i) = 1; conjg(v(1:n-k+i-1)) is stored on exit in A(m-k+i,1:n-k+i-1), and tau in TAU(i).
Definition at line 124 of file zgerq2.f.
subroutine zgerqf | ( | integer | M, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( * ) | TAU, | ||
complex*16, dimension( * ) | WORK, | ||
integer | LWORK, | ||
integer | INFO | ||
) |
ZGERQF
Download ZGERQF + dependencies [TGZ] [ZIP] [TXT]ZGERQF computes an RQ factorization of a complex M-by-N matrix A: A = R * Q.
[in] | M | M is INTEGER The number of rows of the matrix A. M >= 0. |
[in] | N | N is INTEGER The number of columns of the matrix A. N >= 0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the M-by-N matrix A. On exit, if m <= n, the upper triangle of the subarray A(1:m,n-m+1:n) contains the M-by-M upper triangular matrix R; if m >= n, the elements on and above the (m-n)-th subdiagonal contain the M-by-N upper trapezoidal matrix R; the remaining elements, with the array TAU, represent the unitary matrix Q as a product of min(m,n) elementary reflectors (see Further Details). |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | TAU | TAU is COMPLEX*16 array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details). |
[out] | WORK | WORK is COMPLEX*16 array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK(1) returns the optimal LWORK. |
[in] | LWORK | LWORK is INTEGER The dimension of the array WORK. LWORK >= max(1,M). For optimum performance LWORK >= M*NB, where NB is the optimal blocksize. 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. |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value |
The matrix Q is represented as a product of elementary reflectors Q = H(1)**H H(2)**H . . . H(k)**H, where k = min(m,n). Each H(i) has the form H(i) = I - tau * v * v**H where tau is a complex scalar, and v is a complex vector with v(n-k+i+1:n) = 0 and v(n-k+i) = 1; conjg(v(1:n-k+i-1)) is stored on exit in A(m-k+i,1:n-k+i-1), and tau in TAU(i).
Definition at line 139 of file zgerqf.f.
subroutine zgetf2 | ( | integer | M, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
integer, dimension( * ) | IPIV, | ||
integer | INFO | ||
) |
ZGETF2 computes the LU factorization of a general m-by-n matrix using partial pivoting with row interchanges (unblocked algorithm).
Download ZGETF2 + dependencies [TGZ] [ZIP] [TXT]ZGETF2 computes an LU factorization of a general m-by-n matrix A using partial pivoting with row interchanges. The factorization has the form A = P * L * U where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n). This is the right-looking Level 2 BLAS version of the algorithm.
[in] | M | M is INTEGER The number of rows of the matrix A. M >= 0. |
[in] | N | N is INTEGER The number of columns of the matrix A. N >= 0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the m by n matrix to be factored. On exit, the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | IPIV | IPIV is INTEGER array, dimension (min(M,N)) The pivot indices; for 1 <= i <= min(M,N), row i of the matrix was interchanged with row IPIV(i). |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -k, the k-th argument had an illegal value > 0: if INFO = k, U(k,k) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations. |
Definition at line 109 of file zgetf2.f.
subroutine zgetrf | ( | integer | M, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
integer, dimension( * ) | IPIV, | ||
integer | INFO | ||
) |
ZGETRF
Download ZGETRF + dependencies [TGZ] [ZIP] [TXT]ZGETRF computes an LU factorization of a general M-by-N matrix A using partial pivoting with row interchanges. The factorization has the form A = P * L * U where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n). This is the right-looking Level 3 BLAS version of the algorithm.
[in] | M | M is INTEGER The number of rows of the matrix A. M >= 0. |
[in] | N | N is INTEGER The number of columns of the matrix A. N >= 0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the M-by-N matrix to be factored. On exit, the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | IPIV | IPIV is INTEGER array, dimension (min(M,N)) The pivot indices; for 1 <= i <= min(M,N), row i of the matrix was interchanged with row IPIV(i). |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value > 0: if INFO = i, U(i,i) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations. |
Definition at line 109 of file zgetrf.f.
subroutine zgetri | ( | integer | N, |
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
integer, dimension( * ) | IPIV, | ||
complex*16, dimension( * ) | WORK, | ||
integer | LWORK, | ||
integer | INFO | ||
) |
ZGETRI
Download ZGETRI + dependencies [TGZ] [ZIP] [TXT]ZGETRI computes the inverse of a matrix using the LU factorization computed by ZGETRF. This method inverts U and then computes inv(A) by solving the system inv(A)*L = inv(U) for inv(A).
[in] | N | N is INTEGER The order of the matrix A. N >= 0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the factors L and U from the factorization A = P*L*U as computed by ZGETRF. On exit, if INFO = 0, the inverse of the original matrix A. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in] | IPIV | IPIV is INTEGER array, dimension (N) The pivot indices from ZGETRF; for 1<=i<=N, row i of the matrix was interchanged with row IPIV(i). |
[out] | WORK | WORK is COMPLEX*16 array, dimension (MAX(1,LWORK)) On exit, if INFO=0, then WORK(1) returns the optimal LWORK. |
[in] | LWORK | LWORK is INTEGER The dimension of the array WORK. LWORK >= max(1,N). For optimal performance LWORK >= N*NB, where NB is the optimal blocksize returned by ILAENV. 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. |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value > 0: if INFO = i, U(i,i) is exactly zero; the matrix is singular and its inverse could not be computed. |
Definition at line 115 of file zgetri.f.
subroutine zgetrs | ( | character | TRANS, |
integer | N, | ||
integer | NRHS, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
integer, dimension( * ) | IPIV, | ||
complex*16, dimension( ldb, * ) | B, | ||
integer | LDB, | ||
integer | INFO | ||
) |
ZGETRS
Download ZGETRS + dependencies [TGZ] [ZIP] [TXT]ZGETRS solves a system of linear equations A * X = B, A**T * X = B, or A**H * X = B with a general N-by-N matrix A using the LU factorization computed by ZGETRF.
[in] | TRANS | TRANS is 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) |
[in] | N | N is INTEGER The order of the matrix A. N >= 0. |
[in] | NRHS | NRHS is INTEGER The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0. |
[in] | A | A is COMPLEX*16 array, dimension (LDA,N) The factors L and U from the factorization A = P*L*U as computed by ZGETRF. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in] | IPIV | IPIV is INTEGER array, dimension (N) The pivot indices from ZGETRF; for 1<=i<=N, row i of the matrix was interchanged with row IPIV(i). |
[in,out] | B | B is COMPLEX*16 array, dimension (LDB,NRHS) On entry, the right hand side matrix B. On exit, the solution matrix X. |
[in] | LDB | LDB is INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value |
Definition at line 122 of file zgetrs.f.
subroutine zhgeqz | ( | character | JOB, |
character | COMPQ, | ||
character | COMPZ, | ||
integer | N, | ||
integer | ILO, | ||
integer | IHI, | ||
complex*16, dimension( ldh, * ) | H, | ||
integer | LDH, | ||
complex*16, dimension( ldt, * ) | T, | ||
integer | LDT, | ||
complex*16, dimension( * ) | ALPHA, | ||
complex*16, dimension( * ) | BETA, | ||
complex*16, dimension( ldq, * ) | Q, | ||
integer | LDQ, | ||
complex*16, dimension( ldz, * ) | Z, | ||
integer | LDZ, | ||
complex*16, dimension( * ) | WORK, | ||
integer | LWORK, | ||
double precision, dimension( * ) | RWORK, | ||
integer | INFO | ||
) |
ZHGEQZ
Download ZHGEQZ + dependencies [TGZ] [ZIP] [TXT]ZHGEQZ computes the eigenvalues of a complex matrix pair (H,T), where H is an upper Hessenberg matrix and T is upper triangular, using the single-shift QZ method. Matrix pairs of this type are produced by the reduction to generalized upper Hessenberg form of a complex matrix pair (A,B): A = Q1*H*Z1**H, B = Q1*T*Z1**H, as computed by ZGGHRD. If JOB='S', then the Hessenberg-triangular pair (H,T) is also reduced to generalized Schur form, H = Q*S*Z**H, T = Q*P*Z**H, where Q and Z are unitary matrices and S and P are upper triangular. Optionally, the unitary matrix Q from the generalized Schur factorization may be postmultiplied into an input matrix Q1, and the unitary matrix Z may be postmultiplied into an input matrix Z1. If Q1 and Z1 are the unitary matrices from ZGGHRD that reduced the matrix pair (A,B) to generalized Hessenberg form, then the output matrices Q1*Q and Z1*Z are the unitary factors from the generalized Schur factorization of (A,B): A = (Q1*Q)*S*(Z1*Z)**H, B = (Q1*Q)*P*(Z1*Z)**H. To avoid overflow, eigenvalues of the matrix pair (H,T) (equivalently, of (A,B)) are computed as a pair of complex values (alpha,beta). If beta is nonzero, lambda = alpha / beta is an eigenvalue of the generalized nonsymmetric eigenvalue problem (GNEP) A*x = lambda*B*x and if alpha is nonzero, mu = beta / alpha is an eigenvalue of the alternate form of the GNEP mu*A*y = B*y. The values of alpha and beta for the i-th eigenvalue can be read directly from the generalized Schur form: alpha = S(i,i), beta = P(i,i). Ref: C.B. Moler & G.W. Stewart, "An Algorithm for Generalized Matrix Eigenvalue Problems", SIAM J. Numer. Anal., 10(1973), pp. 241--256.
[in] | JOB | JOB is CHARACTER*1 = 'E': Compute eigenvalues only; = 'S': Computer eigenvalues and the Schur form. |
[in] | COMPQ | COMPQ is CHARACTER*1 = 'N': Left Schur vectors (Q) are not computed; = 'I': Q is initialized to the unit matrix and the matrix Q of left Schur vectors of (H,T) is returned; = 'V': Q must contain a unitary matrix Q1 on entry and the product Q1*Q is returned. |
[in] | COMPZ | COMPZ is CHARACTER*1 = 'N': Right Schur vectors (Z) are not computed; = 'I': Q is initialized to the unit matrix and the matrix Z of right Schur vectors of (H,T) is returned; = 'V': Z must contain a unitary matrix Z1 on entry and the product Z1*Z is returned. |
[in] | N | N is INTEGER The order of the matrices H, T, Q, and Z. N >= 0. |
[in] | ILO | ILO is INTEGER |
[in] | IHI | IHI is INTEGER ILO and IHI mark the rows and columns of H which are in Hessenberg form. It is assumed that A is already upper triangular in rows and columns 1:ILO-1 and IHI+1:N. If N > 0, 1 <= ILO <= IHI <= N; if N = 0, ILO=1 and IHI=0. |
[in,out] | H | H is COMPLEX*16 array, dimension (LDH, N) On entry, the N-by-N upper Hessenberg matrix H. On exit, if JOB = 'S', H contains the upper triangular matrix S from the generalized Schur factorization. If JOB = 'E', the diagonal of H matches that of S, but the rest of H is unspecified. |
[in] | LDH | LDH is INTEGER The leading dimension of the array H. LDH >= max( 1, N ). |
[in,out] | T | T is COMPLEX*16 array, dimension (LDT, N) On entry, the N-by-N upper triangular matrix T. On exit, if JOB = 'S', T contains the upper triangular matrix P from the generalized Schur factorization. If JOB = 'E', the diagonal of T matches that of P, but the rest of T is unspecified. |
[in] | LDT | LDT is INTEGER The leading dimension of the array T. LDT >= max( 1, N ). |
[out] | ALPHA | ALPHA is COMPLEX*16 array, dimension (N) The complex scalars alpha that define the eigenvalues of GNEP. ALPHA(i) = S(i,i) in the generalized Schur factorization. |
[out] | BETA | BETA is COMPLEX*16 array, dimension (N) The real non-negative scalars beta that define the eigenvalues of GNEP. BETA(i) = P(i,i) in the generalized Schur factorization. Together, the quantities alpha = ALPHA(j) and beta = BETA(j) represent the j-th eigenvalue of the matrix pair (A,B), in one of the forms lambda = alpha/beta or mu = beta/alpha. Since either lambda or mu may overflow, they should not, in general, be computed. |
[in,out] | Q | Q is COMPLEX*16 array, dimension (LDQ, N) On entry, if COMPZ = 'V', the unitary matrix Q1 used in the reduction of (A,B) to generalized Hessenberg form. On exit, if COMPZ = 'I', the unitary matrix of left Schur vectors of (H,T), and if COMPZ = 'V', the unitary matrix of left Schur vectors of (A,B). Not referenced if COMPZ = 'N'. |
[in] | LDQ | LDQ is INTEGER The leading dimension of the array Q. LDQ >= 1. If COMPQ='V' or 'I', then LDQ >= N. |
[in,out] | Z | Z is COMPLEX*16 array, dimension (LDZ, N) On entry, if COMPZ = 'V', the unitary matrix Z1 used in the reduction of (A,B) to generalized Hessenberg form. On exit, if COMPZ = 'I', the unitary matrix of right Schur vectors of (H,T), and if COMPZ = 'V', the unitary matrix of right Schur vectors of (A,B). Not referenced if COMPZ = 'N'. |
[in] | LDZ | LDZ is INTEGER The leading dimension of the array Z. LDZ >= 1. If COMPZ='V' or 'I', then LDZ >= N. |
[out] | WORK | WORK is COMPLEX*16 array, dimension (MAX(1,LWORK)) On exit, if INFO >= 0, WORK(1) returns the optimal LWORK. |
[in] | LWORK | LWORK is INTEGER The dimension of the array WORK. LWORK >= max(1,N). 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. |
[out] | RWORK | RWORK is DOUBLE PRECISION array, dimension (N) |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value = 1,...,N: the QZ iteration did not converge. (H,T) is not in Schur form, but ALPHA(i) and BETA(i), i=INFO+1,...,N should be correct. = N+1,...,2*N: the shift calculation failed. (H,T) is not in Schur form, but ALPHA(i) and BETA(i), i=INFO-N+1,...,N should be correct. |
We assume that complex ABS works as long as its value is less than overflow.
Definition at line 283 of file zhgeqz.f.
subroutine zla_geamv | ( | integer | TRANS, |
integer | M, | ||
integer | N, | ||
double precision | ALPHA, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( * ) | X, | ||
integer | INCX, | ||
double precision | BETA, | ||
double precision, dimension( * ) | Y, | ||
integer | INCY | ||
) |
ZLA_GEAMV computes a matrix-vector product using a general matrix to calculate error bounds.
Download ZLA_GEAMV + dependencies [TGZ] [ZIP] [TXT]ZLA_GEAMV performs one of the matrix-vector operations y := alpha*abs(A)*abs(x) + beta*abs(y), or y := alpha*abs(A)**T*abs(x) + beta*abs(y), where alpha and beta are scalars, x and y are vectors and A is an m by n matrix. This function is primarily used in calculating error bounds. To protect against underflow during evaluation, components in the resulting vector are perturbed away from zero by (N+1) times the underflow threshold. To prevent unnecessarily large errors for block-structure embedded in general matrices, "symbolically" zero components are not perturbed. A zero entry is considered "symbolic" if all multiplications involved in computing that entry have at least one zero multiplicand.
[in] | TRANS | TRANS is INTEGER On entry, TRANS specifies the operation to be performed as follows: BLAS_NO_TRANS y := alpha*abs(A)*abs(x) + beta*abs(y) BLAS_TRANS y := alpha*abs(A**T)*abs(x) + beta*abs(y) BLAS_CONJ_TRANS y := alpha*abs(A**T)*abs(x) + beta*abs(y) Unchanged on exit. |
[in] | M | M is INTEGER On entry, M specifies the number of rows of the matrix A. M must be at least zero. Unchanged on exit. |
[in] | N | N is INTEGER On entry, N specifies the number of columns of the matrix A. N must be at least zero. Unchanged on exit. |
[in] | ALPHA | ALPHA is DOUBLE PRECISION On entry, ALPHA specifies the scalar alpha. Unchanged on exit. |
[in] | A | A is COMPLEX*16 array of DIMENSION ( LDA, n ) Before entry, the leading m by n part of the array A must contain the matrix of coefficients. Unchanged on exit. |
[in] | LDA | LDA is INTEGER On entry, LDA specifies the first dimension of A as declared in the calling (sub) program. LDA must be at least max( 1, m ). Unchanged on exit. |
[in] | X | X is COMPLEX*16 array of DIMENSION at least ( 1 + ( n - 1 )*abs( INCX ) ) when TRANS = 'N' or 'n' and at least ( 1 + ( m - 1 )*abs( INCX ) ) otherwise. Before entry, the incremented array X must contain the vector x. Unchanged on exit. |
[in] | INCX | INCX is INTEGER On entry, INCX specifies the increment for the elements of X. INCX must not be zero. Unchanged on exit. |
[in] | BETA | BETA is DOUBLE PRECISION On entry, BETA specifies the scalar beta. When BETA is supplied as zero then Y need not be set on input. Unchanged on exit. |
[in,out] | Y | Y is DOUBLE PRECISION array, dimension ( 1 + ( m - 1 )*abs( INCY ) ) when TRANS = 'N' or 'n' and at least ( 1 + ( n - 1 )*abs( INCY ) ) otherwise. Before entry with BETA non-zero, the incremented array Y must contain the vector y. On exit, Y is overwritten by the updated vector y. |
[in] | INCY | INCY is INTEGER On entry, INCY specifies the increment for the elements of Y. INCY must not be zero. Unchanged on exit. Level 2 Blas routine. |
Definition at line 175 of file zla_geamv.f.
DOUBLE PRECISION function zla_gercond_c | ( | character | TRANS, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( ldaf, * ) | AF, | ||
integer | LDAF, | ||
integer, dimension( * ) | IPIV, | ||
double precision, dimension( * ) | C, | ||
logical | CAPPLY, | ||
integer | INFO, | ||
complex*16, dimension( * ) | WORK, | ||
double precision, dimension( * ) | RWORK | ||
) |
ZLA_GERCOND_C computes the infinity norm condition number of op(A)*inv(diag(c)) for general matrices.
Download ZLA_GERCOND_C + dependencies [TGZ] [ZIP] [TXT]ZLA_GERCOND_C computes the infinity norm condition number of op(A) * inv(diag(C)) where C is a DOUBLE PRECISION vector.
[in] | TRANS | TRANS is 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) |
[in] | N | N is INTEGER The number of linear equations, i.e., the order of the matrix A. N >= 0. |
[in] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the N-by-N matrix A |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in] | AF | AF is COMPLEX*16 array, dimension (LDAF,N) The factors L and U from the factorization A = P*L*U as computed by ZGETRF. |
[in] | LDAF | LDAF is INTEGER The leading dimension of the array AF. LDAF >= max(1,N). |
[in] | IPIV | IPIV is INTEGER array, dimension (N) The pivot indices from the factorization A = P*L*U as computed by ZGETRF; row i of the matrix was interchanged with row IPIV(i). |
[in] | C | C is DOUBLE PRECISION array, dimension (N) The vector C in the formula op(A) * inv(diag(C)). |
[in] | CAPPLY | CAPPLY is LOGICAL If .TRUE. then access the vector C in the formula above. |
[out] | INFO | INFO is INTEGER = 0: Successful exit. i > 0: The ith argument is invalid. |
[in] | WORK | WORK is COMPLEX*16 array, dimension (2*N). Workspace. |
[in] | RWORK | RWORK is DOUBLE PRECISION array, dimension (N). Workspace. |
Definition at line 142 of file zla_gercond_c.f.
DOUBLE PRECISION function zla_gercond_x | ( | character | TRANS, |
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( ldaf, * ) | AF, | ||
integer | LDAF, | ||
integer, dimension( * ) | IPIV, | ||
complex*16, dimension( * ) | X, | ||
integer | INFO, | ||
complex*16, dimension( * ) | WORK, | ||
double precision, dimension( * ) | RWORK | ||
) |
ZLA_GERCOND_X computes the infinity norm condition number of op(A)*diag(x) for general matrices.
Download ZLA_GERCOND_X + dependencies [TGZ] [ZIP] [TXT]ZLA_GERCOND_X computes the infinity norm condition number of op(A) * diag(X) where X is a COMPLEX*16 vector.
[in] | TRANS | TRANS is 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) |
[in] | N | N is INTEGER The number of linear equations, i.e., the order of the matrix A. N >= 0. |
[in] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the N-by-N matrix A. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in] | AF | AF is COMPLEX*16 array, dimension (LDAF,N) The factors L and U from the factorization A = P*L*U as computed by ZGETRF. |
[in] | LDAF | LDAF is INTEGER The leading dimension of the array AF. LDAF >= max(1,N). |
[in] | IPIV | IPIV is INTEGER array, dimension (N) The pivot indices from the factorization A = P*L*U as computed by ZGETRF; row i of the matrix was interchanged with row IPIV(i). |
[in] | X | X is COMPLEX*16 array, dimension (N) The vector X in the formula op(A) * diag(X). |
[out] | INFO | INFO is INTEGER = 0: Successful exit. i > 0: The ith argument is invalid. |
[in] | WORK | WORK is COMPLEX*16 array, dimension (2*N). Workspace. |
[in] | RWORK | RWORK is DOUBLE PRECISION array, dimension (N). Workspace. |
Definition at line 135 of file zla_gercond_x.f.
subroutine zla_gerfsx_extended | ( | integer | PREC_TYPE, |
integer | TRANS_TYPE, | ||
integer | N, | ||
integer | NRHS, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( ldaf, * ) | AF, | ||
integer | LDAF, | ||
integer, dimension( * ) | IPIV, | ||
logical | COLEQU, | ||
double precision, dimension( * ) | C, | ||
complex*16, dimension( ldb, * ) | B, | ||
integer | LDB, | ||
complex*16, dimension( ldy, * ) | Y, | ||
integer | LDY, | ||
double precision, dimension( * ) | BERR_OUT, | ||
integer | N_NORMS, | ||
double precision, dimension( nrhs, * ) | ERRS_N, | ||
double precision, dimension( nrhs, * ) | ERRS_C, | ||
complex*16, dimension( * ) | RES, | ||
double precision, dimension( * ) | AYB, | ||
complex*16, dimension( * ) | DY, | ||
complex*16, dimension( * ) | Y_TAIL, | ||
double precision | RCOND, | ||
integer | ITHRESH, | ||
double precision | RTHRESH, | ||
double precision | DZ_UB, | ||
logical | IGNORE_CWISE, | ||
integer | INFO | ||
) |
ZLA_GERFSX_EXTENDED
Download ZLA_GERFSX_EXTENDED + dependencies [TGZ] [ZIP] [TXT]ZLA_GERFSX_EXTENDED improves the computed solution to a system of linear equations by performing extra-precise iterative refinement and provides error bounds and backward error estimates for the solution. This subroutine is called by ZGERFSX to perform iterative refinement. In addition to normwise error bound, the code provides maximum componentwise error bound if possible. See comments for ERRS_N and ERRS_C for details of the error bounds. Note that this subroutine is only resonsible for setting the second fields of ERRS_N and ERRS_C.
[in] | PREC_TYPE | PREC_TYPE is INTEGER Specifies the intermediate precision to be used in refinement. The value is defined by ILAPREC(P) where P is a CHARACTER and P = 'S': Single = 'D': Double = 'I': Indigenous = 'X', 'E': Extra |
[in] | TRANS_TYPE | TRANS_TYPE is INTEGER Specifies the transposition operation on A. The value is defined by ILATRANS(T) where T is a CHARACTER and T = 'N': No transpose = 'T': Transpose = 'C': Conjugate transpose |
[in] | N | N is INTEGER The number of linear equations, i.e., the order of the matrix A. N >= 0. |
[in] | NRHS | NRHS is INTEGER The number of right-hand-sides, i.e., the number of columns of the matrix B. |
[in] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the N-by-N matrix A. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in] | AF | AF is COMPLEX*16 array, dimension (LDAF,N) The factors L and U from the factorization A = P*L*U as computed by ZGETRF. |
[in] | LDAF | LDAF is INTEGER The leading dimension of the array AF. LDAF >= max(1,N). |
[in] | IPIV | IPIV is INTEGER array, dimension (N) The pivot indices from the factorization A = P*L*U as computed by ZGETRF; row i of the matrix was interchanged with row IPIV(i). |
[in] | COLEQU | COLEQU is LOGICAL If .TRUE. then column equilibration was done to A before calling this routine. This is needed to compute the solution and error bounds correctly. |
[in] | C | C is DOUBLE PRECISION array, dimension (N) The column scale factors for A. If COLEQU = .FALSE., C is not accessed. If C is input, each element of C should be a power of the radix to ensure a reliable solution and error estimates. Scaling by powers of the radix does not cause rounding errors unless the result underflows or overflows. Rounding errors during scaling lead to refining with a matrix that is not equivalent to the input matrix, producing error estimates that may not be reliable. |
[in] | B | B is COMPLEX*16 array, dimension (LDB,NRHS) The right-hand-side matrix B. |
[in] | LDB | LDB is INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in,out] | Y | Y is COMPLEX*16 array, dimension (LDY,NRHS) On entry, the solution matrix X, as computed by ZGETRS. On exit, the improved solution matrix Y. |
[in] | LDY | LDY is INTEGER The leading dimension of the array Y. LDY >= max(1,N). |
[out] | BERR_OUT | BERR_OUT is DOUBLE PRECISION array, dimension (NRHS) On exit, BERR_OUT(j) contains the componentwise relative backward error for right-hand-side j from the formula max(i) ( abs(RES(i)) / ( abs(op(A_s))*abs(Y) + abs(B_s) )(i) ) where abs(Z) is the componentwise absolute value of the matrix or vector Z. This is computed by ZLA_LIN_BERR. |
[in] | N_NORMS | N_NORMS is INTEGER Determines which error bounds to return (see ERRS_N and ERRS_C). If N_NORMS >= 1 return normwise error bounds. If N_NORMS >= 2 return componentwise error bounds. |
[in,out] | ERRS_N | ERRS_N is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS) For each right-hand side, this array contains information about various error bounds and condition numbers corresponding to the normwise relative error, which is defined as follows: Normwise relative error in the ith solution vector: max_j (abs(XTRUE(j,i) - X(j,i))) ------------------------------ max_j abs(X(j,i)) The array is indexed by the type of error information as described below. There currently are up to three pieces of information returned. The first index in ERRS_N(i,:) corresponds to the ith right-hand side. The second index in ERRS_N(:,err) contains the following three fields: err = 1 "Trust/don't trust" boolean. Trust the answer if the reciprocal condition number is less than the threshold sqrt(n) * slamch('Epsilon'). err = 2 "Guaranteed" error bound: The estimated forward error, almost certainly within a factor of 10 of the true error so long as the next entry is greater than the threshold sqrt(n) * slamch('Epsilon'). This error bound should only be trusted if the previous boolean is true. err = 3 Reciprocal condition number: Estimated normwise reciprocal condition number. Compared with the threshold sqrt(n) * slamch('Epsilon') to determine if the error estimate is "guaranteed". These reciprocal condition numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some appropriately scaled matrix Z. Let Z = S*A, where S scales each row by a power of the radix so all absolute row sums of Z are approximately 1. This subroutine is only responsible for setting the second field above. See Lapack Working Note 165 for further details and extra cautions. |
[in,out] | ERRS_C | ERRS_C is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS) For each right-hand side, this array contains information about various error bounds and condition numbers corresponding to the componentwise relative error, which is defined as follows: Componentwise relative error in the ith solution vector: abs(XTRUE(j,i) - X(j,i)) max_j ---------------------- abs(X(j,i)) The array is indexed by the right-hand side i (on which the componentwise relative error depends), and the type of error information as described below. There currently are up to three pieces of information returned for each right-hand side. If componentwise accuracy is not requested (PARAMS(3) = 0.0), then ERRS_C is not accessed. If N_ERR_BNDS .LT. 3, then at most the first (:,N_ERR_BNDS) entries are returned. The first index in ERRS_C(i,:) corresponds to the ith right-hand side. The second index in ERRS_C(:,err) contains the following three fields: err = 1 "Trust/don't trust" boolean. Trust the answer if the reciprocal condition number is less than the threshold sqrt(n) * slamch('Epsilon'). err = 2 "Guaranteed" error bound: The estimated forward error, almost certainly within a factor of 10 of the true error so long as the next entry is greater than the threshold sqrt(n) * slamch('Epsilon'). This error bound should only be trusted if the previous boolean is true. err = 3 Reciprocal condition number: Estimated componentwise reciprocal condition number. Compared with the threshold sqrt(n) * slamch('Epsilon') to determine if the error estimate is "guaranteed". These reciprocal condition numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some appropriately scaled matrix Z. Let Z = S*(A*diag(x)), where x is the solution for the current right-hand side and S scales each row of A*diag(x) by a power of the radix so all absolute row sums of Z are approximately 1. This subroutine is only responsible for setting the second field above. See Lapack Working Note 165 for further details and extra cautions. |
[in] | RES | RES is COMPLEX*16 array, dimension (N) Workspace to hold the intermediate residual. |
[in] | AYB | AYB is DOUBLE PRECISION array, dimension (N) Workspace. |
[in] | DY | DY is COMPLEX*16 array, dimension (N) Workspace to hold the intermediate solution. |
[in] | Y_TAIL | Y_TAIL is COMPLEX*16 array, dimension (N) Workspace to hold the trailing bits of the intermediate solution. |
[in] | RCOND | RCOND is DOUBLE PRECISION Reciprocal scaled condition number. This is an estimate of the reciprocal Skeel condition number of the matrix A after equilibration (if done). If this is less than the machine precision (in particular, if it is zero), the matrix is singular to working precision. Note that the error may still be small even if this number is very small and the matrix appears ill- conditioned. |
[in] | ITHRESH | ITHRESH is INTEGER The maximum number of residual computations allowed for refinement. The default is 10. For 'aggressive' set to 100 to permit convergence using approximate factorizations or factorizations other than LU. If the factorization uses a technique other than Gaussian elimination, the guarantees in ERRS_N and ERRS_C may no longer be trustworthy. |
[in] | RTHRESH | RTHRESH is DOUBLE PRECISION Determines when to stop refinement if the error estimate stops decreasing. Refinement will stop when the next solution no longer satisfies norm(dx_{i+1}) < RTHRESH * norm(dx_i) where norm(Z) is the infinity norm of Z. RTHRESH satisfies 0 < RTHRESH <= 1. The default value is 0.5. For 'aggressive' set to 0.9 to permit convergence on extremely ill-conditioned matrices. See LAWN 165 for more details. |
[in] | DZ_UB | DZ_UB is DOUBLE PRECISION Determines when to start considering componentwise convergence. Componentwise convergence is only considered after each component of the solution Y is stable, which we definte as the relative change in each component being less than DZ_UB. The default value is 0.25, requiring the first bit to be stable. See LAWN 165 for more details. |
[in] | IGNORE_CWISE | IGNORE_CWISE is LOGICAL If .TRUE. then ignore componentwise convergence. Default value is .FALSE.. |
[out] | INFO | INFO is INTEGER = 0: Successful exit. < 0: if INFO = -i, the ith argument to ZGETRS had an illegal value |
Definition at line 394 of file zla_gerfsx_extended.f.
DOUBLE PRECISION function zla_gerpvgrw | ( | integer | N, |
integer | NCOLS, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( ldaf, * ) | AF, | ||
integer | LDAF | ||
) |
ZLA_GERPVGRW multiplies a square real matrix by a complex matrix.
Download ZLA_GERPVGRW + dependencies [TGZ] [ZIP] [TXT]ZLA_GERPVGRW computes the reciprocal pivot growth factor norm(A)/norm(U). The "max absolute element" norm is used. If this is much less than 1, the stability of the LU factorization of the (equilibrated) matrix A could be poor. This also means that the solution X, estimated condition numbers, and error bounds could be unreliable.
[in] | N | N is INTEGER The number of linear equations, i.e., the order of the matrix A. N >= 0. |
[in] | NCOLS | NCOLS is INTEGER The number of columns of the matrix A. NCOLS >= 0. |
[in] | A | A is DOUBLE PRECISION array, dimension (LDA,N) On entry, the N-by-N matrix A. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in] | AF | AF is DOUBLE PRECISION array, dimension (LDAF,N) The factors L and U from the factorization A = P*L*U as computed by ZGETRF. |
[in] | LDAF | LDAF is INTEGER The leading dimension of the array AF. LDAF >= max(1,N). |
Definition at line 100 of file zla_gerpvgrw.f.
subroutine ztgevc | ( | character | SIDE, |
character | HOWMNY, | ||
logical, dimension( * ) | SELECT, | ||
integer | N, | ||
complex*16, dimension( lds, * ) | S, | ||
integer | LDS, | ||
complex*16, dimension( ldp, * ) | P, | ||
integer | LDP, | ||
complex*16, dimension( ldvl, * ) | VL, | ||
integer | LDVL, | ||
complex*16, dimension( ldvr, * ) | VR, | ||
integer | LDVR, | ||
integer | MM, | ||
integer | M, | ||
complex*16, dimension( * ) | WORK, | ||
double precision, dimension( * ) | RWORK, | ||
integer | INFO | ||
) |
ZTGEVC
Download ZTGEVC + dependencies [TGZ] [ZIP] [TXT]ZTGEVC computes some or all of the right and/or left eigenvectors of a pair of complex matrices (S,P), where S and P are upper triangular. Matrix pairs of this type are produced by the generalized Schur factorization of a complex matrix pair (A,B): A = Q*S*Z**H, B = Q*P*Z**H as computed by ZGGHRD + ZHGEQZ. The right eigenvector x and the left eigenvector y of (S,P) corresponding to an eigenvalue w are defined by: S*x = w*P*x, (y**H)*S = w*(y**H)*P, where y**H denotes the conjugate tranpose of y. The eigenvalues are not input to this routine, but are computed directly from the diagonal elements of S and P. This routine returns the matrices X and/or Y of right and left eigenvectors of (S,P), or the products Z*X and/or Q*Y, where Z and Q are input matrices. If Q and Z are the unitary factors from the generalized Schur factorization of a matrix pair (A,B), then Z*X and Q*Y are the matrices of right and left eigenvectors of (A,B).
[in] | SIDE | SIDE is CHARACTER*1 = 'R': compute right eigenvectors only; = 'L': compute left eigenvectors only; = 'B': compute both right and left eigenvectors. |
[in] | HOWMNY | HOWMNY is 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, specified by the logical array SELECT. |
[in] | SELECT | SELECT is LOGICAL array, dimension (N) If HOWMNY='S', SELECT specifies the eigenvectors to be computed. The eigenvector corresponding to the j-th eigenvalue is computed if SELECT(j) = .TRUE.. Not referenced if HOWMNY = 'A' or 'B'. |
[in] | N | N is INTEGER The order of the matrices S and P. N >= 0. |
[in] | S | S is COMPLEX*16 array, dimension (LDS,N) The upper triangular matrix S from a generalized Schur factorization, as computed by ZHGEQZ. |
[in] | LDS | LDS is INTEGER The leading dimension of array S. LDS >= max(1,N). |
[in] | P | P is COMPLEX*16 array, dimension (LDP,N) The upper triangular matrix P from a generalized Schur factorization, as computed by ZHGEQZ. P must have real diagonal elements. |
[in] | LDP | LDP is INTEGER The leading dimension of array P. LDP >= max(1,N). |
[in,out] | VL | VL is COMPLEX*16 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 unitary matrix Q of left Schur vectors returned by ZHGEQZ). On exit, if SIDE = 'L' or 'B', VL contains: if HOWMNY = 'A', the matrix Y of left eigenvectors of (S,P); if HOWMNY = 'B', the matrix Q*Y; if HOWMNY = 'S', the left eigenvectors of (S,P) specified by SELECT, stored consecutively in the columns of VL, in the same order as their eigenvalues. Not referenced if SIDE = 'R'. |
[in] | LDVL | LDVL is INTEGER The leading dimension of array VL. LDVL >= 1, and if SIDE = 'L' or 'l' or 'B' or 'b', LDVL >= N. |
[in,out] | VR | VR is COMPLEX*16 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 unitary matrix Z of right Schur vectors returned by ZHGEQZ). On exit, if SIDE = 'R' or 'B', VR contains: if HOWMNY = 'A', the matrix X of right eigenvectors of (S,P); if HOWMNY = 'B', the matrix Z*X; if HOWMNY = 'S', the right eigenvectors of (S,P) specified by SELECT, stored consecutively in the columns of VR, in the same order as their eigenvalues. Not referenced if SIDE = 'L'. |
[in] | LDVR | LDVR is INTEGER The leading dimension of the array VR. LDVR >= 1, and if SIDE = 'R' or 'B', LDVR >= N. |
[in] | MM | MM is INTEGER The number of columns in the arrays VL and/or VR. MM >= M. |
[out] | M | M is 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 eigenvector occupies one column. |
[out] | WORK | WORK is COMPLEX*16 array, dimension (2*N) |
[out] | RWORK | RWORK is DOUBLE PRECISION array, dimension (2*N) |
[out] | INFO | INFO is INTEGER = 0: successful exit. < 0: if INFO = -i, the i-th argument had an illegal value. |
Definition at line 219 of file ztgevc.f.
subroutine ztgexc | ( | logical | WANTQ, |
logical | WANTZ, | ||
integer | N, | ||
complex*16, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex*16, dimension( ldb, * ) | B, | ||
integer | LDB, | ||
complex*16, dimension( ldq, * ) | Q, | ||
integer | LDQ, | ||
complex*16, dimension( ldz, * ) | Z, | ||
integer | LDZ, | ||
integer | IFST, | ||
integer | ILST, | ||
integer | INFO | ||
) |
ZTGEXC
Download ZTGEXC + dependencies [TGZ] [ZIP] [TXT]ZTGEXC reorders the generalized Schur decomposition of a complex matrix pair (A,B), using an unitary equivalence transformation (A, B) := Q * (A, B) * Z**H, so that the diagonal block of (A, B) with row index IFST is moved to row ILST. (A, B) must be in generalized Schur canonical form, that is, A and B are both upper triangular. Optionally, the matrices Q and Z of generalized Schur vectors are updated. Q(in) * A(in) * Z(in)**H = Q(out) * A(out) * Z(out)**H Q(in) * B(in) * Z(in)**H = Q(out) * B(out) * Z(out)**H
[in] | WANTQ | WANTQ is LOGICAL .TRUE. : update the left transformation matrix Q; .FALSE.: do not update Q. |
[in] | WANTZ | WANTZ is LOGICAL .TRUE. : update the right transformation matrix Z; .FALSE.: do not update Z. |
[in] | N | N is INTEGER The order of the matrices A and B. N >= 0. |
[in,out] | A | A is COMPLEX*16 array, dimension (LDA,N) On entry, the upper triangular matrix A in the pair (A, B). On exit, the updated matrix A. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | B | B is COMPLEX*16 array, dimension (LDB,N) On entry, the upper triangular matrix B in the pair (A, B). On exit, the updated matrix B. |
[in] | LDB | LDB is INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in,out] | Q | Q is COMPLEX*16 array, dimension (LDZ,N) On entry, if WANTQ = .TRUE., the unitary matrix Q. On exit, the updated matrix Q. If WANTQ = .FALSE., Q is not referenced. |
[in] | LDQ | LDQ is INTEGER The leading dimension of the array Q. LDQ >= 1; If WANTQ = .TRUE., LDQ >= N. |
[in,out] | Z | Z is COMPLEX*16 array, dimension (LDZ,N) On entry, if WANTZ = .TRUE., the unitary matrix Z. On exit, the updated matrix Z. If WANTZ = .FALSE., Z is not referenced. |
[in] | LDZ | LDZ is INTEGER The leading dimension of the array Z. LDZ >= 1; If WANTZ = .TRUE., LDZ >= N. |
[in] | IFST | IFST is INTEGER |
[in,out] | ILST | ILST is INTEGER Specify the reordering of the diagonal blocks of (A, B). The block with row index IFST is moved to row ILST, by a sequence of swapping between adjacent blocks. |
[out] | INFO | INFO is INTEGER =0: Successful exit. <0: if INFO = -i, the i-th argument had an illegal value. =1: The transformed matrix pair (A, B) would be too far from generalized Schur form; the problem is ill- conditioned. (A, B) may have been partially reordered, and ILST points to the first row of the current position of the block being moved. |
Definition at line 200 of file ztgexc.f.