3 * =========== DOCUMENTATION ===========
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
9 *> Download DGGQRF + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dggqrf.f">
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dggqrf.f">
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dggqrf.f">
21 * SUBROUTINE DGGQRF( N, M, P, A, LDA, TAUA, B, LDB, TAUB, WORK,
24 * .. Scalar Arguments ..
25 * INTEGER INFO, LDA, LDB, LWORK, M, N, P
27 * .. Array Arguments ..
28 * DOUBLE PRECISION A( LDA, * ), B( LDB, * ), TAUA( * ), TAUB( * ),
38 *> DGGQRF computes a generalized QR factorization of an N-by-M matrix A
39 *> and an N-by-P matrix B:
41 *> A = Q*R, B = Q*T*Z,
43 *> where Q is an N-by-N orthogonal matrix, Z is a P-by-P orthogonal
44 *> matrix, and R and T assume one of the forms:
46 *> if N >= M, R = ( R11 ) M , or if N < M, R = ( R11 R12 ) N,
50 *> where R11 is upper triangular, and
52 *> if N <= P, T = ( 0 T12 ) N, or if N > P, T = ( T11 ) N-P,
56 *> where T12 or T21 is upper triangular.
58 *> In particular, if B is square and nonsingular, the GQR factorization
59 *> of A and B implicitly gives the QR factorization of inv(B)*A:
61 *> inv(B)*A = Z**T*(inv(T)*R)
63 *> where inv(B) denotes the inverse of the matrix B, and Z**T denotes the
64 *> transpose of the matrix Z.
73 *> The number of rows of the matrices A and B. N >= 0.
79 *> The number of columns of the matrix A. M >= 0.
85 *> The number of columns of the matrix B. P >= 0.
90 *> A is DOUBLE PRECISION array, dimension (LDA,M)
91 *> On entry, the N-by-M matrix A.
92 *> On exit, the elements on and above the diagonal of the array
93 *> contain the min(N,M)-by-M upper trapezoidal matrix R (R is
94 *> upper triangular if N >= M); the elements below the diagonal,
95 *> with the array TAUA, represent the orthogonal matrix Q as a
96 *> product of min(N,M) elementary reflectors (see Further
103 *> The leading dimension of the array A. LDA >= max(1,N).
108 *> TAUA is DOUBLE PRECISION array, dimension (min(N,M))
109 *> The scalar factors of the elementary reflectors which
110 *> represent the orthogonal matrix Q (see Further Details).
115 *> B is DOUBLE PRECISION array, dimension (LDB,P)
116 *> On entry, the N-by-P matrix B.
117 *> On exit, if N <= P, the upper triangle of the subarray
118 *> B(1:N,P-N+1:P) contains the N-by-N upper triangular matrix T;
119 *> if N > P, the elements on and above the (N-P)-th subdiagonal
120 *> contain the N-by-P upper trapezoidal matrix T; the remaining
121 *> elements, with the array TAUB, represent the orthogonal
122 *> matrix Z as a product of elementary reflectors (see Further
129 *> The leading dimension of the array B. LDB >= max(1,N).
134 *> TAUB is DOUBLE PRECISION array, dimension (min(N,P))
135 *> The scalar factors of the elementary reflectors which
136 *> represent the orthogonal matrix Z (see Further Details).
141 *> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
142 *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
148 *> The dimension of the array WORK. LWORK >= max(1,N,M,P).
149 *> For optimum performance LWORK >= max(N,M,P)*max(NB1,NB2,NB3),
150 *> where NB1 is the optimal blocksize for the QR factorization
151 *> of an N-by-M matrix, NB2 is the optimal blocksize for the
152 *> RQ factorization of an N-by-P matrix, and NB3 is the optimal
153 *> blocksize for a call of DORMQR.
155 *> If LWORK = -1, then a workspace query is assumed; the routine
156 *> only calculates the optimal size of the WORK array, returns
157 *> this value as the first entry of the WORK array, and no error
158 *> message related to LWORK is issued by XERBLA.
164 *> = 0: successful exit
165 *> < 0: if INFO = -i, the i-th argument had an illegal value.
171 *> \author Univ. of Tennessee
172 *> \author Univ. of California Berkeley
173 *> \author Univ. of Colorado Denver
176 *> \date November 2011
178 *> \ingroup doubleOTHERcomputational
180 *> \par Further Details:
181 * =====================
185 *> The matrix Q is represented as a product of elementary reflectors
187 *> Q = H(1) H(2) . . . H(k), where k = min(n,m).
189 *> Each H(i) has the form
191 *> H(i) = I - taua * v * v**T
193 *> where taua is a real scalar, and v is a real vector with
194 *> v(1:i-1) = 0 and v(i) = 1; v(i+1:n) is stored on exit in A(i+1:n,i),
195 *> and taua in TAUA(i).
196 *> To form Q explicitly, use LAPACK subroutine DORGQR.
197 *> To use Q to update another matrix, use LAPACK subroutine DORMQR.
199 *> The matrix Z is represented as a product of elementary reflectors
201 *> Z = H(1) H(2) . . . H(k), where k = min(n,p).
203 *> Each H(i) has the form
205 *> H(i) = I - taub * v * v**T
207 *> where taub is a real scalar, and v is a real vector with
208 *> v(p-k+i+1:p) = 0 and v(p-k+i) = 1; v(1:p-k+i-1) is stored on exit in
209 *> B(n-k+i,1:p-k+i-1), and taub in TAUB(i).
210 *> To form Z explicitly, use LAPACK subroutine DORGRQ.
211 *> To use Z to update another matrix, use LAPACK subroutine DORMRQ.
214 * =====================================================================
215 SUBROUTINE DGGQRF( N, M, P, A, LDA, TAUA, B, LDB, TAUB, WORK,
218 * -- LAPACK computational routine (version 3.4.0) --
219 * -- LAPACK is a software package provided by Univ. of Tennessee, --
220 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
223 * .. Scalar Arguments ..
224 INTEGER INFO, LDA, LDB, LWORK, M, N, P
226 * .. Array Arguments ..
227 DOUBLE PRECISION A( LDA, * ), B( LDB, * ), TAUA( * ), TAUB( * ),
231 * =====================================================================
233 * .. Local Scalars ..
235 INTEGER LOPT, LWKOPT, NB, NB1, NB2, NB3
237 * .. External Subroutines ..
238 EXTERNAL DGEQRF, DGERQF, DORMQR, XERBLA
240 * .. External Functions ..
244 * .. Intrinsic Functions ..
245 INTRINSIC INT, MAX, MIN
247 * .. Executable Statements ..
249 * Test the input parameters
252 NB1 = ILAENV( 1, 'DGEQRF', ' ', N, M, -1, -1 )
253 NB2 = ILAENV( 1, 'DGERQF', ' ', N, P, -1, -1 )
254 NB3 = ILAENV( 1, 'DORMQR', ' ', N, M, P, -1 )
255 NB = MAX( NB1, NB2, NB3 )
256 LWKOPT = MAX( N, M, P )*NB
258 LQUERY = ( LWORK.EQ.-1 )
261 ELSE IF( M.LT.0 ) THEN
263 ELSE IF( P.LT.0 ) THEN
265 ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
267 ELSE IF( LDB.LT.MAX( 1, N ) ) THEN
269 ELSE IF( LWORK.LT.MAX( 1, N, M, P ) .AND. .NOT.LQUERY ) THEN
273 CALL XERBLA( 'DGGQRF', -INFO )
275 ELSE IF( LQUERY ) THEN
279 * QR factorization of N-by-M matrix A: A = Q*R
281 CALL DGEQRF( N, M, A, LDA, TAUA, WORK, LWORK, INFO )
284 * Update B := Q**T*B.
286 CALL DORMQR( 'Left', 'Transpose', N, P, MIN( N, M ), A, LDA, TAUA,
287 $ B, LDB, WORK, LWORK, INFO )
288 LOPT = MAX( LOPT, INT( WORK( 1 ) ) )
290 * RQ factorization of N-by-P matrix B: B = T*Z.
292 CALL DGERQF( N, P, B, LDB, TAUB, WORK, LWORK, INFO )
293 WORK( 1 ) = MAX( LOPT, INT( WORK( 1 ) ) )