ggqrf - Man Page

ggqrf: Generalized QR factor

Synopsis

Functions

subroutine cggqrf (n, m, p, a, lda, taua, b, ldb, taub, work, lwork, info)
CGGQRF
subroutine dggqrf (n, m, p, a, lda, taua, b, ldb, taub, work, lwork, info)
DGGQRF
subroutine sggqrf (n, m, p, a, lda, taua, b, ldb, taub, work, lwork, info)
SGGQRF
subroutine zggqrf (n, m, p, a, lda, taua, b, ldb, taub, work, lwork, info)
ZGGQRF

Detailed Description

Function Documentation

subroutine cggqrf (integer n, integer m, integer p, complex, dimension( lda, * ) a, integer lda, complex, dimension( * ) taua, complex, dimension( ldb, * ) b, integer ldb, complex, dimension( * ) taub, complex, dimension( * ) work, integer lwork, integer info)

CGGQRF  

Purpose:

 CGGQRF computes a generalized QR factorization of an N-by-M matrix A
 and an N-by-P matrix B:

             A = Q*R,        B = Q*T*Z,

 where Q is an N-by-N unitary matrix, Z is a P-by-P unitary matrix,
 and R and T assume one of the forms:

 if N >= M,  R = ( R11 ) M  ,   or if N < M,  R = ( R11  R12 ) N,
                 (  0  ) N-M                         N   M-N
                    M

 where R11 is upper triangular, and

 if N <= P,  T = ( 0  T12 ) N,   or if N > P,  T = ( T11 ) N-P,
                  P-N  N                           ( T21 ) P
                                                      P

 where T12 or T21 is upper triangular.

 In particular, if B is square and nonsingular, the GQR factorization
 of A and B implicitly gives the QR factorization of inv(B)*A:

              inv(B)*A = Z**H * (inv(T)*R)

 where inv(B) denotes the inverse of the matrix B, and Z' denotes the
 conjugate transpose of matrix Z.
Parameters

N

          N is INTEGER
          The number of rows of the matrices A and B. N >= 0.

M

          M is INTEGER
          The number of columns of the matrix A.  M >= 0.

P

          P is INTEGER
          The number of columns of the matrix B.  P >= 0.

A

          A is COMPLEX array, dimension (LDA,M)
          On entry, the N-by-M matrix A.
          On exit, the elements on and above the diagonal of the array
          contain the min(N,M)-by-M upper trapezoidal matrix R (R is
          upper triangular if N >= M); the elements below the diagonal,
          with the array TAUA, represent the unitary matrix Q as a
          product of min(N,M) elementary reflectors (see Further
          Details).

LDA

          LDA is INTEGER
          The leading dimension of the array A. LDA >= max(1,N).

TAUA

          TAUA is COMPLEX array, dimension (min(N,M))
          The scalar factors of the elementary reflectors which
          represent the unitary matrix Q (see Further Details).

B

          B is COMPLEX array, dimension (LDB,P)
          On entry, the N-by-P matrix B.
          On exit, if N <= P, the upper triangle of the subarray
          B(1:N,P-N+1:P) contains the N-by-N upper triangular matrix T;
          if N > P, the elements on and above the (N-P)-th subdiagonal
          contain the N-by-P upper trapezoidal matrix T; the remaining
          elements, with the array TAUB, represent the unitary
          matrix Z as a product of elementary reflectors (see Further
          Details).

LDB

          LDB is INTEGER
          The leading dimension of the array B. LDB >= max(1,N).

TAUB

          TAUB is COMPLEX array, dimension (min(N,P))
          The scalar factors of the elementary reflectors which
          represent the unitary matrix Z (see Further Details).

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,M,P).
          For optimum performance LWORK >= max(N,M,P)*max(NB1,NB2,NB3),
          where NB1 is the optimal blocksize for the QR factorization
          of an N-by-M matrix, NB2 is the optimal blocksize for the
          RQ factorization of an N-by-P matrix, and NB3 is the optimal
          blocksize for a call of CUNMQR.

          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

Univ. of Tennessee

Univ. of California Berkeley

Univ. of Colorado Denver

NAG Ltd.

Further Details:

  The matrix Q is represented as a product of elementary reflectors

     Q = H(1) H(2) . . . H(k), where k = min(n,m).

  Each H(i) has the form

     H(i) = I - taua * v * v**H

  where taua is a complex scalar, and v is a complex vector with
  v(1:i-1) = 0 and v(i) = 1; v(i+1:n) is stored on exit in A(i+1:n,i),
  and taua in TAUA(i).
  To form Q explicitly, use LAPACK subroutine CUNGQR.
  To use Q to update another matrix, use LAPACK subroutine CUNMQR.

  The matrix Z is represented as a product of elementary reflectors

     Z = H(1) H(2) . . . H(k), where k = min(n,p).

  Each H(i) has the form

     H(i) = I - taub * v * v**H

  where taub is a complex scalar, and v is a complex vector with
  v(p-k+i+1:p) = 0 and v(p-k+i) = 1; v(1:p-k+i-1) is stored on exit in
  B(n-k+i,1:p-k+i-1), and taub in TAUB(i).
  To form Z explicitly, use LAPACK subroutine CUNGRQ.
  To use Z to update another matrix, use LAPACK subroutine CUNMRQ.

Definition at line 213 of file cggqrf.f.

subroutine dggqrf (integer n, integer m, integer p, double precision, dimension( lda, * ) a, integer lda, double precision, dimension( * ) taua, double precision, dimension( ldb, * ) b, integer ldb, double precision, dimension( * ) taub, double precision, dimension( * ) work, integer lwork, integer info)

DGGQRF  

Purpose:

 DGGQRF computes a generalized QR factorization of an N-by-M matrix A
 and an N-by-P matrix B:

             A = Q*R,        B = Q*T*Z,

 where Q is an N-by-N orthogonal matrix, Z is a P-by-P orthogonal
 matrix, and R and T assume one of the forms:

 if N >= M,  R = ( R11 ) M  ,   or if N < M,  R = ( R11  R12 ) N,
                 (  0  ) N-M                         N   M-N
                    M

 where R11 is upper triangular, and

 if N <= P,  T = ( 0  T12 ) N,   or if N > P,  T = ( T11 ) N-P,
                  P-N  N                           ( T21 ) P
                                                      P

 where T12 or T21 is upper triangular.

 In particular, if B is square and nonsingular, the GQR factorization
 of A and B implicitly gives the QR factorization of inv(B)*A:

              inv(B)*A = Z**T*(inv(T)*R)

 where inv(B) denotes the inverse of the matrix B, and Z**T denotes the
 transpose of the matrix Z.
Parameters

N

          N is INTEGER
          The number of rows of the matrices A and B. N >= 0.

M

          M is INTEGER
          The number of columns of the matrix A.  M >= 0.

P

          P is INTEGER
          The number of columns of the matrix B.  P >= 0.

A

          A is DOUBLE PRECISION array, dimension (LDA,M)
          On entry, the N-by-M matrix A.
          On exit, the elements on and above the diagonal of the array
          contain the min(N,M)-by-M upper trapezoidal matrix R (R is
          upper triangular if N >= M); the elements below the diagonal,
          with the array TAUA, represent the orthogonal matrix Q as a
          product of min(N,M) elementary reflectors (see Further
          Details).

LDA

          LDA is INTEGER
          The leading dimension of the array A. LDA >= max(1,N).

TAUA

          TAUA is DOUBLE PRECISION array, dimension (min(N,M))
          The scalar factors of the elementary reflectors which
          represent the orthogonal matrix Q (see Further Details).

B

          B is DOUBLE PRECISION array, dimension (LDB,P)
          On entry, the N-by-P matrix B.
          On exit, if N <= P, the upper triangle of the subarray
          B(1:N,P-N+1:P) contains the N-by-N upper triangular matrix T;
          if N > P, the elements on and above the (N-P)-th subdiagonal
          contain the N-by-P upper trapezoidal matrix T; the remaining
          elements, with the array TAUB, represent the orthogonal
          matrix Z as a product of elementary reflectors (see Further
          Details).

LDB

          LDB is INTEGER
          The leading dimension of the array B. LDB >= max(1,N).

TAUB

          TAUB is DOUBLE PRECISION array, dimension (min(N,P))
          The scalar factors of the elementary reflectors which
          represent the orthogonal matrix Z (see Further Details).

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,M,P).
          For optimum performance LWORK >= max(N,M,P)*max(NB1,NB2,NB3),
          where NB1 is the optimal blocksize for the QR factorization
          of an N-by-M matrix, NB2 is the optimal blocksize for the
          RQ factorization of an N-by-P matrix, and NB3 is the optimal
          blocksize for a call of DORMQR.

          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

Univ. of Tennessee

Univ. of California Berkeley

Univ. of Colorado Denver

NAG Ltd.

Further Details:

  The matrix Q is represented as a product of elementary reflectors

     Q = H(1) H(2) . . . H(k), where k = min(n,m).

  Each H(i) has the form

     H(i) = I - taua * v * v**T

  where taua is a real scalar, and v is a real vector with
  v(1:i-1) = 0 and v(i) = 1; v(i+1:n) is stored on exit in A(i+1:n,i),
  and taua in TAUA(i).
  To form Q explicitly, use LAPACK subroutine DORGQR.
  To use Q to update another matrix, use LAPACK subroutine DORMQR.

  The matrix Z is represented as a product of elementary reflectors

     Z = H(1) H(2) . . . H(k), where k = min(n,p).

  Each H(i) has the form

     H(i) = I - taub * v * v**T

  where taub is a real scalar, and v is a real vector with
  v(p-k+i+1:p) = 0 and v(p-k+i) = 1; v(1:p-k+i-1) is stored on exit in
  B(n-k+i,1:p-k+i-1), and taub in TAUB(i).
  To form Z explicitly, use LAPACK subroutine DORGRQ.
  To use Z to update another matrix, use LAPACK subroutine DORMRQ.

Definition at line 213 of file dggqrf.f.

subroutine sggqrf (integer n, integer m, integer p, real, dimension( lda, * ) a, integer lda, real, dimension( * ) taua, real, dimension( ldb, * ) b, integer ldb, real, dimension( * ) taub, real, dimension( * ) work, integer lwork, integer info)

SGGQRF  

Purpose:

 SGGQRF computes a generalized QR factorization of an N-by-M matrix A
 and an N-by-P matrix B:

             A = Q*R,        B = Q*T*Z,

 where Q is an N-by-N orthogonal matrix, Z is a P-by-P orthogonal
 matrix, and R and T assume one of the forms:

 if N >= M,  R = ( R11 ) M  ,   or if N < M,  R = ( R11  R12 ) N,
                 (  0  ) N-M                         N   M-N
                    M

 where R11 is upper triangular, and

 if N <= P,  T = ( 0  T12 ) N,   or if N > P,  T = ( T11 ) N-P,
                  P-N  N                           ( T21 ) P
                                                      P

 where T12 or T21 is upper triangular.

 In particular, if B is square and nonsingular, the GQR factorization
 of A and B implicitly gives the QR factorization of inv(B)*A:

              inv(B)*A = Z**T*(inv(T)*R)

 where inv(B) denotes the inverse of the matrix B, and Z**T denotes the
 transpose of the matrix Z.
Parameters

N

          N is INTEGER
          The number of rows of the matrices A and B. N >= 0.

M

          M is INTEGER
          The number of columns of the matrix A.  M >= 0.

P

          P is INTEGER
          The number of columns of the matrix B.  P >= 0.

A

          A is REAL array, dimension (LDA,M)
          On entry, the N-by-M matrix A.
          On exit, the elements on and above the diagonal of the array
          contain the min(N,M)-by-M upper trapezoidal matrix R (R is
          upper triangular if N >= M); the elements below the diagonal,
          with the array TAUA, represent the orthogonal matrix Q as a
          product of min(N,M) elementary reflectors (see Further
          Details).

LDA

          LDA is INTEGER
          The leading dimension of the array A. LDA >= max(1,N).

TAUA

          TAUA is REAL array, dimension (min(N,M))
          The scalar factors of the elementary reflectors which
          represent the orthogonal matrix Q (see Further Details).

B

          B is REAL array, dimension (LDB,P)
          On entry, the N-by-P matrix B.
          On exit, if N <= P, the upper triangle of the subarray
          B(1:N,P-N+1:P) contains the N-by-N upper triangular matrix T;
          if N > P, the elements on and above the (N-P)-th subdiagonal
          contain the N-by-P upper trapezoidal matrix T; the remaining
          elements, with the array TAUB, represent the orthogonal
          matrix Z as a product of elementary reflectors (see Further
          Details).

LDB

          LDB is INTEGER
          The leading dimension of the array B. LDB >= max(1,N).

TAUB

          TAUB is REAL array, dimension (min(N,P))
          The scalar factors of the elementary reflectors which
          represent the orthogonal matrix Z (see Further Details).

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,M,P).
          For optimum performance LWORK >= max(N,M,P)*max(NB1,NB2,NB3),
          where NB1 is the optimal blocksize for the QR factorization
          of an N-by-M matrix, NB2 is the optimal blocksize for the
          RQ factorization of an N-by-P matrix, and NB3 is the optimal
          blocksize for a call of SORMQR.

          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

Univ. of Tennessee

Univ. of California Berkeley

Univ. of Colorado Denver

NAG Ltd.

Further Details:

  The matrix Q is represented as a product of elementary reflectors

     Q = H(1) H(2) . . . H(k), where k = min(n,m).

  Each H(i) has the form

     H(i) = I - taua * v * v**T

  where taua is a real scalar, and v is a real vector with
  v(1:i-1) = 0 and v(i) = 1; v(i+1:n) is stored on exit in A(i+1:n,i),
  and taua in TAUA(i).
  To form Q explicitly, use LAPACK subroutine SORGQR.
  To use Q to update another matrix, use LAPACK subroutine SORMQR.

  The matrix Z is represented as a product of elementary reflectors

     Z = H(1) H(2) . . . H(k), where k = min(n,p).

  Each H(i) has the form

     H(i) = I - taub * v * v**T

  where taub is a real scalar, and v is a real vector with
  v(p-k+i+1:p) = 0 and v(p-k+i) = 1; v(1:p-k+i-1) is stored on exit in
  B(n-k+i,1:p-k+i-1), and taub in TAUB(i).
  To form Z explicitly, use LAPACK subroutine SORGRQ.
  To use Z to update another matrix, use LAPACK subroutine SORMRQ.

Definition at line 213 of file sggqrf.f.

subroutine zggqrf (integer n, integer m, integer p, complex*16, dimension( lda, * ) a, integer lda, complex*16, dimension( * ) taua, complex*16, dimension( ldb, * ) b, integer ldb, complex*16, dimension( * ) taub, complex*16, dimension( * ) work, integer lwork, integer info)

ZGGQRF  

Purpose:

 ZGGQRF computes a generalized QR factorization of an N-by-M matrix A
 and an N-by-P matrix B:

             A = Q*R,        B = Q*T*Z,

 where Q is an N-by-N unitary matrix, Z is a P-by-P unitary matrix,
 and R and T assume one of the forms:

 if N >= M,  R = ( R11 ) M  ,   or if N < M,  R = ( R11  R12 ) N,
                 (  0  ) N-M                         N   M-N
                    M

 where R11 is upper triangular, and

 if N <= P,  T = ( 0  T12 ) N,   or if N > P,  T = ( T11 ) N-P,
                  P-N  N                           ( T21 ) P
                                                      P

 where T12 or T21 is upper triangular.

 In particular, if B is square and nonsingular, the GQR factorization
 of A and B implicitly gives the QR factorization of inv(B)*A:

              inv(B)*A = Z**H * (inv(T)*R)

 where inv(B) denotes the inverse of the matrix B, and Z**H denotes the
 conjugate transpose of matrix Z.
Parameters

N

          N is INTEGER
          The number of rows of the matrices A and B. N >= 0.

M

          M is INTEGER
          The number of columns of the matrix A.  M >= 0.

P

          P is INTEGER
          The number of columns of the matrix B.  P >= 0.

A

          A is COMPLEX*16 array, dimension (LDA,M)
          On entry, the N-by-M matrix A.
          On exit, the elements on and above the diagonal of the array
          contain the min(N,M)-by-M upper trapezoidal matrix R (R is
          upper triangular if N >= M); the elements below the diagonal,
          with the array TAUA, represent the unitary matrix Q as a
          product of min(N,M) elementary reflectors (see Further
          Details).

LDA

          LDA is INTEGER
          The leading dimension of the array A. LDA >= max(1,N).

TAUA

          TAUA is COMPLEX*16 array, dimension (min(N,M))
          The scalar factors of the elementary reflectors which
          represent the unitary matrix Q (see Further Details).

B

          B is COMPLEX*16 array, dimension (LDB,P)
          On entry, the N-by-P matrix B.
          On exit, if N <= P, the upper triangle of the subarray
          B(1:N,P-N+1:P) contains the N-by-N upper triangular matrix T;
          if N > P, the elements on and above the (N-P)-th subdiagonal
          contain the N-by-P upper trapezoidal matrix T; the remaining
          elements, with the array TAUB, represent the unitary
          matrix Z as a product of elementary reflectors (see Further
          Details).

LDB

          LDB is INTEGER
          The leading dimension of the array B. LDB >= max(1,N).

TAUB

          TAUB is COMPLEX*16 array, dimension (min(N,P))
          The scalar factors of the elementary reflectors which
          represent the unitary matrix Z (see Further Details).

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,M,P).
          For optimum performance LWORK >= max(N,M,P)*max(NB1,NB2,NB3),
          where NB1 is the optimal blocksize for the QR factorization
          of an N-by-M matrix, NB2 is the optimal blocksize for the
          RQ factorization of an N-by-P matrix, and NB3 is the optimal
          blocksize for a call of ZUNMQR.

          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

Univ. of Tennessee

Univ. of California Berkeley

Univ. of Colorado Denver

NAG Ltd.

Further Details:

  The matrix Q is represented as a product of elementary reflectors

     Q = H(1) H(2) . . . H(k), where k = min(n,m).

  Each H(i) has the form

     H(i) = I - taua * v * v**H

  where taua is a complex scalar, and v is a complex vector with
  v(1:i-1) = 0 and v(i) = 1; v(i+1:n) is stored on exit in A(i+1:n,i),
  and taua in TAUA(i).
  To form Q explicitly, use LAPACK subroutine ZUNGQR.
  To use Q to update another matrix, use LAPACK subroutine ZUNMQR.

  The matrix Z is represented as a product of elementary reflectors

     Z = H(1) H(2) . . . H(k), where k = min(n,p).

  Each H(i) has the form

     H(i) = I - taub * v * v**H

  where taub is a complex scalar, and v is a complex vector with
  v(p-k+i+1:p) = 0 and v(p-k+i) = 1; v(1:p-k+i-1) is stored on exit in
  B(n-k+i,1:p-k+i-1), and taub in TAUB(i).
  To form Z explicitly, use LAPACK subroutine ZUNGRQ.
  To use Z to update another matrix, use LAPACK subroutine ZUNMRQ.

Definition at line 213 of file zggqrf.f.

Author

Generated automatically by Doxygen for LAPACK from the source code.

Info

Tue Nov 28 2023 12:08:43 Version 3.12.0 LAPACK