Provided by: liblapack-doc_3.12.0-3build1.1_all 
      
    
NAME
       laqz3 - laqz3: step in ggev3, gges3
SYNOPSIS
   Functions
       subroutine claqz3 (ilschur, ilq, ilz, n, ilo, ihi, nshifts, nblock_desired, alpha, beta, a, lda, b, ldb,
           q, ldq, z, ldz, qc, ldqc, zc, ldzc, work, lwork, info)
           CLAQZ3
       recursive subroutine dlaqz3 (ilschur, ilq, ilz, n, ilo, ihi, nw, a, lda, b, ldb, q, ldq, z, ldz, ns, nd,
           alphar, alphai, beta, qc, ldqc, zc, ldzc, work, lwork, rec, info)
           DLAQZ3
       recursive subroutine slaqz3 (ilschur, ilq, ilz, n, ilo, ihi, nw, a, lda, b, ldb, q, ldq, z, ldz, ns, nd,
           alphar, alphai, beta, qc, ldqc, zc, ldzc, work, lwork, rec, info)
           SLAQZ3
       subroutine zlaqz3 (ilschur, ilq, ilz, n, ilo, ihi, nshifts, nblock_desired, alpha, beta, a, lda, b, ldb,
           q, ldq, z, ldz, qc, ldqc, zc, ldzc, work, lwork, info)
           ZLAQZ3
Detailed Description
Function Documentation
   subroutine claqz3 (logical, intent(in) ilschur, logical, intent(in) ilq, logical, intent(in) ilz, integer,
       intent(in) n, integer, intent(in) ilo, integer, intent(in) ihi, integer, intent(in) nshifts, integer,
       intent(in) nblock_desired, complex, dimension( * ), intent(inout) alpha, complex, dimension( * ),
       intent(inout) beta, complex, dimension( lda, * ), intent(inout) a, integer, intent(in) lda, complex,
       dimension( ldb, * ), intent(inout) b, integer, intent(in) ldb, complex, dimension( ldq, * ),
       intent(inout) q, integer, intent(in) ldq, complex, dimension( ldz, * ), intent(inout) z, integer,
       intent(in) ldz, complex, dimension( ldqc, * ), intent(inout) qc, integer, intent(in) ldqc, complex,
       dimension( ldzc, * ), intent(inout) zc, integer, intent(in) ldzc, complex, dimension( * ), intent(inout)
       work, integer, intent(in) lwork, integer, intent(out) info)
       CLAQZ3
       Purpose:
            CLAQZ3 Executes a single multishift QZ sweep
       Parameters
           ILSCHUR
                     ILSCHUR is LOGICAL
                         Determines whether or not to update the full Schur form
           ILQ
                     ILQ is LOGICAL
                         Determines whether or not to update the matrix Q
           ILZ
                     ILZ is LOGICAL
                         Determines whether or not to update the matrix Z
           N
                     N is INTEGER
                     The order of the matrices A, B, Q, and Z.  N >= 0.
           ILO
                     ILO is INTEGER
           IHI
                     IHI is INTEGER
           NSHIFTS
                     NSHIFTS is INTEGER
                     The desired number of shifts to use
           NBLOCK_DESIRED
                     NBLOCK_DESIRED is INTEGER
                     The desired size of the computational windows
           ALPHA
                     ALPHA is COMPLEX array. SR contains
                     the alpha parts of the shifts to use.
           BETA
                     BETA is COMPLEX array. SS contains
                     the scale of the shifts to use.
           A
                     A is COMPLEX array, dimension (LDA, N)
           LDA
                     LDA is INTEGER
                     The leading dimension of the array A.  LDA >= max( 1, N ).
           B
                     B is COMPLEX array, dimension (LDB, N)
           LDB
                     LDB is INTEGER
                     The leading dimension of the array B.  LDB >= max( 1, N ).
           Q
                     Q is COMPLEX array, dimension (LDQ, N)
           LDQ
                     LDQ is INTEGER
           Z
                     Z is COMPLEX array, dimension (LDZ, N)
           LDZ
                     LDZ is INTEGER
           QC
                     QC is COMPLEX array, dimension (LDQC, NBLOCK_DESIRED)
           LDQC
                     LDQC is INTEGER
           ZC
                     ZC is COMPLEX array, dimension (LDZC, NBLOCK_DESIRED)
           LDZC
                     LDZ is INTEGER
           WORK
                     WORK is COMPLEX array, dimension (MAX(1,LWORK))
                     On exit, if INFO >= 0, WORK(1) returns the optimal LWORK.
           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.
           INFO
                     INFO is INTEGER
                     = 0: successful exit
                     < 0: if INFO = -i, the i-th argument had an illegal value
       Author
           Thijs Steel, KU Leuven
       Date
           May 2020
   recursive subroutine dlaqz3 (logical, intent(in) ilschur, logical, intent(in) ilq, logical, intent(in) ilz,
       integer, intent(in) n, integer, intent(in) ilo, integer, intent(in) ihi, integer, intent(in) nw, double
       precision, dimension( lda, * ), intent(inout) a, integer, intent(in) lda, double precision, dimension(
       ldb, * ), intent(inout) b, integer, intent(in) ldb, double precision, dimension( ldq, * ), intent(inout)
       q, integer, intent(in) ldq, double precision, dimension( ldz, * ), intent(inout) z, integer, intent(in)
       ldz, integer, intent(out) ns, integer, intent(out) nd, double precision, dimension( * ), intent(inout)
       alphar, double precision, dimension( * ), intent(inout) alphai, double precision, dimension( * ),
       intent(inout) beta, double precision, dimension( ldqc, * ) qc, integer, intent(in) ldqc, double
       precision, dimension( ldzc, * ) zc, integer, intent(in) ldzc, double precision, dimension( * ) work,
       integer, intent(in) lwork, integer, intent(in) rec, integer, intent(out) info)
       DLAQZ3
       Purpose:
            DLAQZ3 performs AED
       Parameters
           ILSCHUR
                     ILSCHUR is LOGICAL
                         Determines whether or not to update the full Schur form
           ILQ
                     ILQ is LOGICAL
                         Determines whether or not to update the matrix Q
           ILZ
                     ILZ is LOGICAL
                         Determines whether or not to update the matrix Z
           N
                     N is INTEGER
                     The order of the matrices A, B, Q, and Z.  N >= 0.
           ILO
                     ILO is INTEGER
           IHI
                     IHI is INTEGER
                     ILO and IHI mark the rows and columns of (A,B) which
                     are to be normalized
           NW
                     NW is INTEGER
                     The desired size of the deflation window.
           A
                     A is DOUBLE PRECISION array, dimension (LDA, N)
           LDA
                     LDA is INTEGER
                     The leading dimension of the array A.  LDA >= max( 1, N ).
           B
                     B is DOUBLE PRECISION array, dimension (LDB, N)
           LDB
                     LDB is INTEGER
                     The leading dimension of the array B.  LDB >= max( 1, N ).
           Q
                     Q is DOUBLE PRECISION array, dimension (LDQ, N)
           LDQ
                     LDQ is INTEGER
           Z
                     Z is DOUBLE PRECISION array, dimension (LDZ, N)
           LDZ
                     LDZ is INTEGER
           NS
                     NS is INTEGER
                     The number of unconverged eigenvalues available to
                     use as shifts.
           ND
                     ND is INTEGER
                     The number of converged eigenvalues found.
           ALPHAR
                     ALPHAR is DOUBLE PRECISION array, dimension (N)
                     The real parts of each scalar alpha defining an eigenvalue
                     of GNEP.
           ALPHAI
                     ALPHAI is DOUBLE PRECISION array, dimension (N)
                     The imaginary parts of each scalar alpha defining an
                     eigenvalue of GNEP.
                     If ALPHAI(j) is zero, then the j-th eigenvalue is real; if
                     positive, then the j-th and (j+1)-st eigenvalues are a
                     complex conjugate pair, with ALPHAI(j+1) = -ALPHAI(j).
           BETA
                     BETA is DOUBLE PRECISION array, dimension (N)
                     The scalars beta that define the eigenvalues of GNEP.
                     Together, the quantities alpha = (ALPHAR(j),ALPHAI(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.
           QC
                     QC is DOUBLE PRECISION array, dimension (LDQC, NW)
           LDQC
                     LDQC is INTEGER
           ZC
                     ZC is DOUBLE PRECISION array, dimension (LDZC, NW)
           LDZC
                     LDZ is INTEGER
           WORK
                     WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
                     On exit, if INFO >= 0, WORK(1) returns the optimal LWORK.
           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.
           REC
                     REC is INTEGER
                        REC indicates the current recursion level. Should be set
                        to 0 on first call.
           INFO
                     INFO is INTEGER
                     = 0: successful exit
                     < 0: if INFO = -i, the i-th argument had an illegal value
       Author
           Thijs Steel, KU Leuven
       Date
           May 2020
   recursive subroutine slaqz3 (logical, intent(in) ilschur, logical, intent(in) ilq, logical, intent(in) ilz,
       integer, intent(in) n, integer, intent(in) ilo, integer, intent(in) ihi, integer, intent(in) nw, real,
       dimension( lda, * ), intent(inout) a, integer, intent(in) lda, real, dimension( ldb, * ), intent(inout)
       b, integer, intent(in) ldb, real, dimension( ldq, * ), intent(inout) q, integer, intent(in) ldq, real,
       dimension( ldz, * ), intent(inout) z, integer, intent(in) ldz, integer, intent(out) ns, integer,
       intent(out) nd, real, dimension( * ), intent(inout) alphar, real, dimension( * ), intent(inout) alphai,
       real, dimension( * ), intent(inout) beta, real, dimension( ldqc, * ) qc, integer, intent(in) ldqc, real,
       dimension( ldzc, * ) zc, integer, intent(in) ldzc, real, dimension( * ) work, integer, intent(in) lwork,
       integer, intent(in) rec, integer, intent(out) info)
       SLAQZ3
       Purpose:
            SLAQZ3 performs AED
       Parameters
           ILSCHUR
                     ILSCHUR is LOGICAL
                         Determines whether or not to update the full Schur form
           ILQ
                     ILQ is LOGICAL
                         Determines whether or not to update the matrix Q
           ILZ
                     ILZ is LOGICAL
                         Determines whether or not to update the matrix Z
           N
                     N is INTEGER
                     The order of the matrices A, B, Q, and Z.  N >= 0.
           ILO
                     ILO is INTEGER
           IHI
                     IHI is INTEGER
                     ILO and IHI mark the rows and columns of (A,B) which
                     are to be normalized
           NW
                     NW is INTEGER
                     The desired size of the deflation window.
           A
                     A is REAL array, dimension (LDA, N)
           LDA
                     LDA is INTEGER
                     The leading dimension of the array A.  LDA >= max( 1, N ).
           B
                     B is REAL array, dimension (LDB, N)
           LDB
                     LDB is INTEGER
                     The leading dimension of the array B.  LDB >= max( 1, N ).
           Q
                     Q is REAL array, dimension (LDQ, N)
           LDQ
                     LDQ is INTEGER
           Z
                     Z is REAL array, dimension (LDZ, N)
           LDZ
                     LDZ is INTEGER
           NS
                     NS is INTEGER
                     The number of unconverged eigenvalues available to
                     use as shifts.
           ND
                     ND is INTEGER
                     The number of converged eigenvalues found.
           ALPHAR
                     ALPHAR is REAL array, dimension (N)
                     The real parts of each scalar alpha defining an eigenvalue
                     of GNEP.
           ALPHAI
                     ALPHAI is REAL array, dimension (N)
                     The imaginary parts of each scalar alpha defining an
                     eigenvalue of GNEP.
                     If ALPHAI(j) is zero, then the j-th eigenvalue is real; if
                     positive, then the j-th and (j+1)-st eigenvalues are a
                     complex conjugate pair, with ALPHAI(j+1) = -ALPHAI(j).
           BETA
                     BETA is REAL array, dimension (N)
                     The scalars beta that define the eigenvalues of GNEP.
                     Together, the quantities alpha = (ALPHAR(j),ALPHAI(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.
           QC
                     QC is REAL array, dimension (LDQC, NW)
           LDQC
                     LDQC is INTEGER
           ZC
                     ZC is REAL array, dimension (LDZC, NW)
           LDZC
                     LDZ is INTEGER
           WORK
                     WORK is REAL array, dimension (MAX(1,LWORK))
                     On exit, if INFO >= 0, WORK(1) returns the optimal LWORK.
           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.
           REC
                     REC is INTEGER
                        REC indicates the current recursion level. Should be set
                        to 0 on first call.
           INFO
                     INFO is INTEGER
                     = 0: successful exit
                     < 0: if INFO = -i, the i-th argument had an illegal value
       Author
           Thijs Steel, KU Leuven
       Date
           May 2020
   subroutine zlaqz3 (logical, intent(in) ilschur, logical, intent(in) ilq, logical, intent(in) ilz, integer,
       intent(in) n, integer, intent(in) ilo, integer, intent(in) ihi, integer, intent(in) nshifts, integer,
       intent(in) nblock_desired, complex*16, dimension( * ), intent(inout) alpha, complex*16, dimension( * ),
       intent(inout) beta, complex*16, dimension( lda, * ), intent(inout) a, integer, intent(in) lda,
       complex*16, dimension( ldb, * ), intent(inout) b, integer, intent(in) ldb, complex*16, dimension( ldq,
       * ), intent(inout) q, integer, intent(in) ldq, complex*16, dimension( ldz, * ), intent(inout) z, integer,
       intent(in) ldz, complex*16, dimension( ldqc, * ), intent(inout) qc, integer, intent(in) ldqc, complex*16,
       dimension( ldzc, * ), intent(inout) zc, integer, intent(in) ldzc, complex*16, dimension( * ),
       intent(inout) work, integer, intent(in) lwork, integer, intent(out) info)
       ZLAQZ3
       Purpose:
            ZLAQZ3 Executes a single multishift QZ sweep
       Parameters
           ILSCHUR
                     ILSCHUR is LOGICAL
                         Determines whether or not to update the full Schur form
           ILQ
                     ILQ is LOGICAL
                         Determines whether or not to update the matrix Q
           ILZ
                     ILZ is LOGICAL
                         Determines whether or not to update the matrix Z
           N
                     N is INTEGER
                     The order of the matrices A, B, Q, and Z.  N >= 0.
           ILO
                     ILO is INTEGER
           IHI
                     IHI is INTEGER
           NSHIFTS
                     NSHIFTS is INTEGER
                     The desired number of shifts to use
           NBLOCK_DESIRED
                     NBLOCK_DESIRED is INTEGER
                     The desired size of the computational windows
           ALPHA
                     ALPHA is COMPLEX*16 array. SR contains
                     the alpha parts of the shifts to use.
           BETA
                     BETA is COMPLEX*16 array. SS contains
                     the scale of the shifts to use.
           A
                     A is COMPLEX*16 array, dimension (LDA, N)
           LDA
                     LDA is INTEGER
                     The leading dimension of the array A.  LDA >= max( 1, N ).
           B
                     B is COMPLEX*16 array, dimension (LDB, N)
           LDB
                     LDB is INTEGER
                     The leading dimension of the array B.  LDB >= max( 1, N ).
           Q
                     Q is COMPLEX*16 array, dimension (LDQ, N)
           LDQ
                     LDQ is INTEGER
           Z
                     Z is COMPLEX*16 array, dimension (LDZ, N)
           LDZ
                     LDZ is INTEGER
           QC
                     QC is COMPLEX*16 array, dimension (LDQC, NBLOCK_DESIRED)
           LDQC
                     LDQC is INTEGER
           ZC
                     ZC is COMPLEX*16 array, dimension (LDZC, NBLOCK_DESIRED)
           LDZC
                     LDZ is INTEGER
           WORK
                     WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
                     On exit, if INFO >= 0, WORK(1) returns the optimal LWORK.
           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.
           INFO
                     INFO is INTEGER
                     = 0: successful exit
                     < 0: if INFO = -i, the i-th argument had an illegal value
       Author
           Thijs Steel, KU Leuven
       Date
           May 2020
Author
       Generated automatically by Doxygen for LAPACK from the source code.
Version 3.12.0                               Fri Aug 9 2024 02:33:22                                    laqz3(3)