LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dgeqr2.f
Go to the documentation of this file.
1*> \brief \b DGEQR2 computes the QR factorization of a general rectangular matrix using an unblocked algorithm.
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download DGEQR2 + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgeqr2.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgeqr2.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgeqr2.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE DGEQR2( M, N, A, LDA, TAU, WORK, INFO )
20*
21* .. Scalar Arguments ..
22* INTEGER INFO, LDA, M, N
23* ..
24* .. Array Arguments ..
25* DOUBLE PRECISION A( LDA, * ), TAU( * ), WORK( * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> DGEQR2 computes a QR factorization of a real m-by-n matrix A:
35*>
36*> A = Q * ( R ),
37*> ( 0 )
38*>
39*> where:
40*>
41*> Q is a m-by-m orthogonal matrix;
42*> R is an upper-triangular n-by-n matrix;
43*> 0 is a (m-n)-by-n zero matrix, if m > n.
44*>
45*> \endverbatim
46*
47* Arguments:
48* ==========
49*
50*> \param[in] M
51*> \verbatim
52*> M is INTEGER
53*> The number of rows of the matrix A. M >= 0.
54*> \endverbatim
55*>
56*> \param[in] N
57*> \verbatim
58*> N is INTEGER
59*> The number of columns of the matrix A. N >= 0.
60*> \endverbatim
61*>
62*> \param[in,out] A
63*> \verbatim
64*> A is DOUBLE PRECISION array, dimension (LDA,N)
65*> On entry, the m by n matrix A.
66*> On exit, the elements on and above the diagonal of the array
67*> contain the min(m,n) by n upper trapezoidal matrix R (R is
68*> upper triangular if m >= n); the elements below the diagonal,
69*> with the array TAU, represent the orthogonal matrix Q as a
70*> product of elementary reflectors (see Further Details).
71*> \endverbatim
72*>
73*> \param[in] LDA
74*> \verbatim
75*> LDA is INTEGER
76*> The leading dimension of the array A. LDA >= max(1,M).
77*> \endverbatim
78*>
79*> \param[out] TAU
80*> \verbatim
81*> TAU is DOUBLE PRECISION array, dimension (min(M,N))
82*> The scalar factors of the elementary reflectors (see Further
83*> Details).
84*> \endverbatim
85*>
86*> \param[out] WORK
87*> \verbatim
88*> WORK is DOUBLE PRECISION array, dimension (N)
89*> \endverbatim
90*>
91*> \param[out] INFO
92*> \verbatim
93*> INFO is INTEGER
94*> = 0: successful exit
95*> < 0: if INFO = -i, the i-th argument had an illegal value
96*> \endverbatim
97*
98* Authors:
99* ========
100*
101*> \author Univ. of Tennessee
102*> \author Univ. of California Berkeley
103*> \author Univ. of Colorado Denver
104*> \author NAG Ltd.
105*
106*> \ingroup geqr2
107*
108*> \par Further Details:
109* =====================
110*>
111*> \verbatim
112*>
113*> The matrix Q is represented as a product of elementary reflectors
114*>
115*> Q = H(1) H(2) . . . H(k), where k = min(m,n).
116*>
117*> Each H(i) has the form
118*>
119*> H(i) = I - tau * v * v**T
120*>
121*> where tau is a real scalar, and v is a real vector with
122*> v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i),
123*> and tau in TAU(i).
124*> \endverbatim
125*>
126* =====================================================================
127 SUBROUTINE dgeqr2( M, N, A, LDA, TAU, WORK, INFO )
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*
192 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgeqr2(m, n, a, lda, tau, work, info)
DGEQR2 computes the QR factorization of a general rectangular matrix using an unblocked algorithm.
Definition dgeqr2.f:128
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