LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ dgeqr2()

subroutine dgeqr2 ( integer m,
integer n,
double precision, dimension( lda, * ) a,
integer lda,
double precision, dimension( * ) tau,
double precision, dimension( * ) work,
integer info )

DGEQR2 computes the QR factorization of a general rectangular matrix using an unblocked algorithm.

Download DGEQR2 + dependencies [TGZ] [ZIP] [TXT]

Purpose:
!>
!> DGEQR2 computes a QR factorization of a real m-by-n matrix A:
!>
!>    A = Q * ( R ),
!>            ( 0 )
!>
!> where:
!>
!>    Q is a m-by-m orthogonal matrix;
!>    R is an upper-triangular n-by-n matrix;
!>    0 is a (m-n)-by-n zero matrix, if m > n.
!>
!> 
Parameters
[in]M
!>          M is INTEGER
!>          The number of rows of the matrix A.  M >= 0.
!> 
[in]N
!>          N is INTEGER
!>          The number of columns of the matrix A.  N >= 0.
!> 
[in,out]A
!>          A is DOUBLE PRECISION array, dimension (LDA,N)
!>          On entry, the m by n matrix A.
!>          On exit, the elements on and above the diagonal of the array
!>          contain the min(m,n) by n upper trapezoidal matrix R (R is
!>          upper triangular if m >= n); the elements below the diagonal,
!>          with the array TAU, represent the orthogonal matrix Q as a
!>          product of elementary reflectors (see Further Details).
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,M).
!> 
[out]TAU
!>          TAU is DOUBLE PRECISION array, dimension (min(M,N))
!>          The scalar factors of the elementary reflectors (see Further
!>          Details).
!> 
[out]WORK
!>          WORK is DOUBLE PRECISION array, dimension (N)
!> 
[out]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(m,n).
!>
!>  Each H(i) has the form
!>
!>     H(i) = I - tau * v * v**T
!>
!>  where tau is a real scalar, and v is a real vector with
!>  v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i),
!>  and tau in TAU(i).
!> 

Definition at line 127 of file dgeqr2.f.

128*
129* -- LAPACK computational routine --
130* -- LAPACK is a software package provided by Univ. of Tennessee, --
131* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
132*
133* .. Scalar Arguments ..
134 INTEGER INFO, LDA, M, N
135* ..
136* .. Array Arguments ..
137 DOUBLE PRECISION A( LDA, * ), TAU( * ), WORK( * )
138* ..
139*
140* =====================================================================
141*
142* .. Parameters ..
143 DOUBLE PRECISION ONE
144 parameter( one = 1.0d+0 )
145* ..
146* .. Local Scalars ..
147 INTEGER I, K
148* ..
149* .. External Subroutines ..
150 EXTERNAL dlarf1f, dlarfg, xerbla
151* ..
152* .. Intrinsic Functions ..
153 INTRINSIC max, min
154* ..
155* .. Executable Statements ..
156*
157* Test the input arguments
158*
159 info = 0
160 IF( m.LT.0 ) THEN
161 info = -1
162 ELSE IF( n.LT.0 ) THEN
163 info = -2
164 ELSE IF( lda.LT.max( 1, m ) ) THEN
165 info = -4
166 END IF
167 IF( info.NE.0 ) THEN
168 CALL xerbla( 'DGEQR2', -info )
169 RETURN
170 END IF
171*
172 k = min( m, n )
173*
174 DO 10 i = 1, k
175*
176* Generate elementary reflector H(i) to annihilate A(i+1:m,i)
177*
178 CALL dlarfg( m-i+1, a( i, i ), a( min( i+1, m ), i ), 1,
179 $ tau( i ) )
180 IF( i.LT.n ) THEN
181*
182* Apply H(i) to A(i:m,i+1:n) from the left
183*
184 CALL dlarf1f( 'Left', m-i+1, n-i, a( i, i ), 1, tau( i ),
185 $ a( i, i+1 ), lda, work )
186 END IF
187 10 CONTINUE
188 RETURN
189*
190* End of DGEQR2
191*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dlarf1f(side, m, n, v, incv, tau, c, ldc, work)
DLARF1F applies an elementary reflector to a general rectangular
Definition dlarf1f.f:157
subroutine dlarfg(n, alpha, x, incx, tau)
DLARFG generates an elementary reflector (Householder matrix).
Definition dlarfg.f:104
Here is the call graph for this function:
Here is the caller graph for this function: