1
SUBROUTINE ZGELQ2( M, N, A, LDA, TAU, WORK, INFO )
3
* -- LAPACK routine (version 3.0) --
4
* Univ. of Tennessee, Univ. of California Berkeley, NAG Ltd.,
5
* Courant Institute, Argonne National Lab, and Rice University
8
* .. Scalar Arguments ..
9
INTEGER INFO, LDA, M, N
11
* .. Array Arguments ..
12
COMPLEX*16 A( LDA, * ), TAU( * ), WORK( * )
18
* ZGELQ2 computes an LQ factorization of a complex m by n matrix A:
25
* The number of rows of the matrix A. M >= 0.
28
* The number of columns of the matrix A. N >= 0.
30
* A (input/output) COMPLEX*16 array, dimension (LDA,N)
31
* On entry, the m by n matrix A.
32
* On exit, the elements on and below the diagonal of the array
33
* contain the m by min(m,n) lower trapezoidal matrix L (L is
34
* lower triangular if m <= n); the elements above the diagonal,
35
* with the array TAU, represent the unitary matrix Q as a
36
* product of elementary reflectors (see Further Details).
39
* The leading dimension of the array A. LDA >= max(1,M).
41
* TAU (output) COMPLEX*16 array, dimension (min(M,N))
42
* The scalar factors of the elementary reflectors (see Further
45
* WORK (workspace) COMPLEX*16 array, dimension (M)
47
* INFO (output) INTEGER
48
* = 0: successful exit
49
* < 0: if INFO = -i, the i-th argument had an illegal value
54
* The matrix Q is represented as a product of elementary reflectors
56
* Q = H(k)' . . . H(2)' H(1)', where k = min(m,n).
58
* Each H(i) has the form
60
* H(i) = I - tau * v * v'
62
* where tau is a complex scalar, and v is a complex vector with
63
* v(1:i-1) = 0 and v(i) = 1; conjg(v(i+1:n)) is stored on exit in
64
* A(i,i+1:n), and tau in TAU(i).
66
* =====================================================================
70
PARAMETER ( ONE = ( 1.0D+0, 0.0D+0 ) )
76
* .. External Subroutines ..
77
EXTERNAL XERBLA, ZLACGV, ZLARF, ZLARFG
79
* .. Intrinsic Functions ..
82
* .. Executable Statements ..
84
* Test the input arguments
89
ELSE IF( N.LT.0 ) THEN
91
ELSE IF( LDA.LT.MAX( 1, M ) ) THEN
95
CALL XERBLA( 'ZGELQ2', -INFO )
103
* Generate elementary reflector H(i) to annihilate A(i,i+1:n)
105
CALL ZLACGV( N-I+1, A( I, I ), LDA )
107
CALL ZLARFG( N-I+1, ALPHA, A( I, MIN( I+1, N ) ), LDA,
111
* Apply H(i) to A(i+1:m,i:n) from the right
114
CALL ZLARF( 'Right', M-I, N-I+1, A( I, I ), LDA, TAU( I ),
115
$ A( I+1, I ), LDA, WORK )
118
CALL ZLACGV( N-I+1, A( I, I ), LDA )