001:       SUBROUTINE DGGRQF( M, P, N, A, LDA, TAUA, B, LDB, TAUB, WORK,
002:      $                   LWORK, INFO )
003: *
004: *  -- LAPACK routine (version 3.2) --
005: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
006: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
007: *     November 2006
008: *
009: *     .. Scalar Arguments ..
010:       INTEGER            INFO, LDA, LDB, LWORK, M, N, P
011: *     ..
012: *     .. Array Arguments ..
013:       DOUBLE PRECISION   A( LDA, * ), B( LDB, * ), TAUA( * ), TAUB( * ),
014:      $                   WORK( * )
015: *     ..
016: *
017: *  Purpose
018: *  =======
019: *
020: *  DGGRQF computes a generalized RQ factorization of an M-by-N matrix A
021: *  and a P-by-N matrix B:
022: *
023: *              A = R*Q,        B = Z*T*Q,
024: *
025: *  where Q is an N-by-N orthogonal matrix, Z is a P-by-P orthogonal
026: *  matrix, and R and T assume one of the forms:
027: *
028: *  if M <= N,  R = ( 0  R12 ) M,   or if M > N,  R = ( R11 ) M-N,
029: *                   N-M  M                           ( R21 ) N
030: *                                                       N
031: *
032: *  where R12 or R21 is upper triangular, and
033: *
034: *  if P >= N,  T = ( T11 ) N  ,   or if P < N,  T = ( T11  T12 ) P,
035: *                  (  0  ) P-N                         P   N-P
036: *                     N
037: *
038: *  where T11 is upper triangular.
039: *
040: *  In particular, if B is square and nonsingular, the GRQ factorization
041: *  of A and B implicitly gives the RQ factorization of A*inv(B):
042: *
043: *               A*inv(B) = (R*inv(T))*Z'
044: *
045: *  where inv(B) denotes the inverse of the matrix B, and Z' denotes the
046: *  transpose of the matrix Z.
047: *
048: *  Arguments
049: *  =========
050: *
051: *  M       (input) INTEGER
052: *          The number of rows of the matrix A.  M >= 0.
053: *
054: *  P       (input) INTEGER
055: *          The number of rows of the matrix B.  P >= 0.
056: *
057: *  N       (input) INTEGER
058: *          The number of columns of the matrices A and B. N >= 0.
059: *
060: *  A       (input/output) DOUBLE PRECISION array, dimension (LDA,N)
061: *          On entry, the M-by-N matrix A.
062: *          On exit, if M <= N, the upper triangle of the subarray
063: *          A(1:M,N-M+1:N) contains the M-by-M upper triangular matrix R;
064: *          if M > N, the elements on and above the (M-N)-th subdiagonal
065: *          contain the M-by-N upper trapezoidal matrix R; the remaining
066: *          elements, with the array TAUA, represent the orthogonal
067: *          matrix Q as a product of elementary reflectors (see Further
068: *          Details).
069: *
070: *  LDA     (input) INTEGER
071: *          The leading dimension of the array A. LDA >= max(1,M).
072: *
073: *  TAUA    (output) DOUBLE PRECISION array, dimension (min(M,N))
074: *          The scalar factors of the elementary reflectors which
075: *          represent the orthogonal matrix Q (see Further Details).
076: *
077: *  B       (input/output) DOUBLE PRECISION array, dimension (LDB,N)
078: *          On entry, the P-by-N matrix B.
079: *          On exit, the elements on and above the diagonal of the array
080: *          contain the min(P,N)-by-N upper trapezoidal matrix T (T is
081: *          upper triangular if P >= N); the elements below the diagonal,
082: *          with the array TAUB, represent the orthogonal matrix Z as a
083: *          product of elementary reflectors (see Further Details).
084: *
085: *  LDB     (input) INTEGER
086: *          The leading dimension of the array B. LDB >= max(1,P).
087: *
088: *  TAUB    (output) DOUBLE PRECISION array, dimension (min(P,N))
089: *          The scalar factors of the elementary reflectors which
090: *          represent the orthogonal matrix Z (see Further Details).
091: *
092: *  WORK    (workspace/output) DOUBLE PRECISION array, dimension (MAX(1,LWORK))
093: *          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
094: *
095: *  LWORK   (input) INTEGER
096: *          The dimension of the array WORK. LWORK >= max(1,N,M,P).
097: *          For optimum performance LWORK >= max(N,M,P)*max(NB1,NB2,NB3),
098: *          where NB1 is the optimal blocksize for the RQ factorization
099: *          of an M-by-N matrix, NB2 is the optimal blocksize for the
100: *          QR factorization of a P-by-N matrix, and NB3 is the optimal
101: *          blocksize for a call of DORMRQ.
102: *
103: *          If LWORK = -1, then a workspace query is assumed; the routine
104: *          only calculates the optimal size of the WORK array, returns
105: *          this value as the first entry of the WORK array, and no error
106: *          message related to LWORK is issued by XERBLA.
107: *
108: *  INFO    (output) INTEGER
109: *          = 0:  successful exit
110: *          < 0:  if INF0= -i, the i-th argument had an illegal value.
111: *
112: *  Further Details
113: *  ===============
114: *
115: *  The matrix Q is represented as a product of elementary reflectors
116: *
117: *     Q = H(1) H(2) . . . H(k), where k = min(m,n).
118: *
119: *  Each H(i) has the form
120: *
121: *     H(i) = I - taua * v * v'
122: *
123: *  where taua is a real scalar, and v is a real vector with
124: *  v(n-k+i+1:n) = 0 and v(n-k+i) = 1; v(1:n-k+i-1) is stored on exit in
125: *  A(m-k+i,1:n-k+i-1), and taua in TAUA(i).
126: *  To form Q explicitly, use LAPACK subroutine DORGRQ.
127: *  To use Q to update another matrix, use LAPACK subroutine DORMRQ.
128: *
129: *  The matrix Z is represented as a product of elementary reflectors
130: *
131: *     Z = H(1) H(2) . . . H(k), where k = min(p,n).
132: *
133: *  Each H(i) has the form
134: *
135: *     H(i) = I - taub * v * v'
136: *
137: *  where taub is a real scalar, and v is a real vector with
138: *  v(1:i-1) = 0 and v(i) = 1; v(i+1:p) is stored on exit in B(i+1:p,i),
139: *  and taub in TAUB(i).
140: *  To form Z explicitly, use LAPACK subroutine DORGQR.
141: *  To use Z to update another matrix, use LAPACK subroutine DORMQR.
142: *
143: *  =====================================================================
144: *
145: *     .. Local Scalars ..
146:       LOGICAL            LQUERY
147:       INTEGER            LOPT, LWKOPT, NB, NB1, NB2, NB3
148: *     ..
149: *     .. External Subroutines ..
150:       EXTERNAL           DGEQRF, DGERQF, DORMRQ, XERBLA
151: *     ..
152: *     .. External Functions ..
153:       INTEGER            ILAENV
154:       EXTERNAL           ILAENV
155: *     ..
156: *     .. Intrinsic Functions ..
157:       INTRINSIC          INT, MAX, MIN
158: *     ..
159: *     .. Executable Statements ..
160: *
161: *     Test the input parameters
162: *
163:       INFO = 0
164:       NB1 = ILAENV( 1, 'DGERQF', ' ', M, N, -1, -1 )
165:       NB2 = ILAENV( 1, 'DGEQRF', ' ', P, N, -1, -1 )
166:       NB3 = ILAENV( 1, 'DORMRQ', ' ', M, N, P, -1 )
167:       NB = MAX( NB1, NB2, NB3 )
168:       LWKOPT = MAX( N, M, P )*NB
169:       WORK( 1 ) = LWKOPT
170:       LQUERY = ( LWORK.EQ.-1 )
171:       IF( M.LT.0 ) THEN
172:          INFO = -1
173:       ELSE IF( P.LT.0 ) THEN
174:          INFO = -2
175:       ELSE IF( N.LT.0 ) THEN
176:          INFO = -3
177:       ELSE IF( LDA.LT.MAX( 1, M ) ) THEN
178:          INFO = -5
179:       ELSE IF( LDB.LT.MAX( 1, P ) ) THEN
180:          INFO = -8
181:       ELSE IF( LWORK.LT.MAX( 1, M, P, N ) .AND. .NOT.LQUERY ) THEN
182:          INFO = -11
183:       END IF
184:       IF( INFO.NE.0 ) THEN
185:          CALL XERBLA( 'DGGRQF', -INFO )
186:          RETURN
187:       ELSE IF( LQUERY ) THEN
188:          RETURN
189:       END IF
190: *
191: *     RQ factorization of M-by-N matrix A: A = R*Q
192: *
193:       CALL DGERQF( M, N, A, LDA, TAUA, WORK, LWORK, INFO )
194:       LOPT = WORK( 1 )
195: *
196: *     Update B := B*Q'
197: *
198:       CALL DORMRQ( 'Right', 'Transpose', P, N, MIN( M, N ),
199:      $             A( MAX( 1, M-N+1 ), 1 ), LDA, TAUA, B, LDB, WORK,
200:      $             LWORK, INFO )
201:       LOPT = MAX( LOPT, INT( WORK( 1 ) ) )
202: *
203: *     QR factorization of P-by-N matrix B: B = Z*T
204: *
205:       CALL DGEQRF( P, N, B, LDB, TAUB, WORK, LWORK, INFO )
206:       WORK( 1 ) = MAX( LOPT, INT( WORK( 1 ) ) )
207: *
208:       RETURN
209: *
210: *     End of DGGRQF
211: *
212:       END
213: