LAPACK  3.9.0
LAPACK: Linear Algebra PACKage
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 *> \htmlonly
9 *> Download DGEQR2 + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgeqr2.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgeqr2.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgeqr2.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE DGEQR2( M, N, A, LDA, TAU, WORK, INFO )
22 *
23 * .. Scalar Arguments ..
24 * INTEGER INFO, LDA, M, N
25 * ..
26 * .. Array Arguments ..
27 * DOUBLE PRECISION A( LDA, * ), TAU( * ), WORK( * )
28 * ..
29 *
30 *
31 *> \par Purpose:
32 * =============
33 *>
34 *> \verbatim
35 *>
36 *> DGEQR2 computes a QR factorization of a real m-by-n matrix A:
37 *>
38 *> A = Q * ( R ),
39 *> ( 0 )
40 *>
41 *> where:
42 *>
43 *> Q is a m-by-m orthogonal matrix;
44 *> R is an upper-triangular n-by-n matrix;
45 *> 0 is a (m-n)-by-n zero matrix, if m > n.
46 *>
47 *> \endverbatim
48 *
49 * Arguments:
50 * ==========
51 *
52 *> \param[in] M
53 *> \verbatim
54 *> M is INTEGER
55 *> The number of rows of the matrix A. M >= 0.
56 *> \endverbatim
57 *>
58 *> \param[in] N
59 *> \verbatim
60 *> N is INTEGER
61 *> The number of columns of the matrix A. N >= 0.
62 *> \endverbatim
63 *>
64 *> \param[in,out] A
65 *> \verbatim
66 *> A is DOUBLE PRECISION array, dimension (LDA,N)
67 *> On entry, the m by n matrix A.
68 *> On exit, the elements on and above the diagonal of the array
69 *> contain the min(m,n) by n upper trapezoidal matrix R (R is
70 *> upper triangular if m >= n); the elements below the diagonal,
71 *> with the array TAU, represent the orthogonal matrix Q as a
72 *> product of elementary reflectors (see Further Details).
73 *> \endverbatim
74 *>
75 *> \param[in] LDA
76 *> \verbatim
77 *> LDA is INTEGER
78 *> The leading dimension of the array A. LDA >= max(1,M).
79 *> \endverbatim
80 *>
81 *> \param[out] TAU
82 *> \verbatim
83 *> TAU is DOUBLE PRECISION array, dimension (min(M,N))
84 *> The scalar factors of the elementary reflectors (see Further
85 *> Details).
86 *> \endverbatim
87 *>
88 *> \param[out] WORK
89 *> \verbatim
90 *> WORK is DOUBLE PRECISION array, dimension (N)
91 *> \endverbatim
92 *>
93 *> \param[out] INFO
94 *> \verbatim
95 *> INFO is INTEGER
96 *> = 0: successful exit
97 *> < 0: if INFO = -i, the i-th argument had an illegal value
98 *> \endverbatim
99 *
100 * Authors:
101 * ========
102 *
103 *> \author Univ. of Tennessee
104 *> \author Univ. of California Berkeley
105 *> \author Univ. of Colorado Denver
106 *> \author NAG Ltd.
107 *
108 *> \date November 2019
109 *
110 *> \ingroup doubleGEcomputational
111 *
112 *> \par Further Details:
113 * =====================
114 *>
115 *> \verbatim
116 *>
117 *> The matrix Q is represented as a product of elementary reflectors
118 *>
119 *> Q = H(1) H(2) . . . H(k), where k = min(m,n).
120 *>
121 *> Each H(i) has the form
122 *>
123 *> H(i) = I - tau * v * v**T
124 *>
125 *> where tau is a real scalar, and v is a real vector with
126 *> v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i),
127 *> and tau in TAU(i).
128 *> \endverbatim
129 *>
130 * =====================================================================
131  SUBROUTINE dgeqr2( M, N, A, LDA, TAU, WORK, INFO )
132 *
133 * -- LAPACK computational routine (version 3.9.0) --
134 * -- LAPACK is a software package provided by Univ. of Tennessee, --
135 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
136 * November 2019
137 *
138 * .. Scalar Arguments ..
139  INTEGER INFO, LDA, M, N
140 * ..
141 * .. Array Arguments ..
142  DOUBLE PRECISION A( LDA, * ), TAU( * ), WORK( * )
143 * ..
144 *
145 * =====================================================================
146 *
147 * .. Parameters ..
148  DOUBLE PRECISION ONE
149  parameter( one = 1.0d+0 )
150 * ..
151 * .. Local Scalars ..
152  INTEGER I, K
153  DOUBLE PRECISION AII
154 * ..
155 * .. External Subroutines ..
156  EXTERNAL dlarf, dlarfg, xerbla
157 * ..
158 * .. Intrinsic Functions ..
159  INTRINSIC max, min
160 * ..
161 * .. Executable Statements ..
162 *
163 * Test the input arguments
164 *
165  info = 0
166  IF( m.LT.0 ) THEN
167  info = -1
168  ELSE IF( n.LT.0 ) THEN
169  info = -2
170  ELSE IF( lda.LT.max( 1, m ) ) THEN
171  info = -4
172  END IF
173  IF( info.NE.0 ) THEN
174  CALL xerbla( 'DGEQR2', -info )
175  RETURN
176  END IF
177 *
178  k = min( m, n )
179 *
180  DO 10 i = 1, k
181 *
182 * Generate elementary reflector H(i) to annihilate A(i+1:m,i)
183 *
184  CALL dlarfg( m-i+1, a( i, i ), a( min( i+1, m ), i ), 1,
185  $ tau( i ) )
186  IF( i.LT.n ) THEN
187 *
188 * Apply H(i) to A(i:m,i+1:n) from the left
189 *
190  aii = a( i, i )
191  a( i, i ) = one
192  CALL dlarf( 'Left', m-i+1, n-i, a( i, i ), 1, tau( i ),
193  $ a( i, i+1 ), lda, work )
194  a( i, i ) = aii
195  END IF
196  10 CONTINUE
197  RETURN
198 *
199 * End of DGEQR2
200 *
201  END
dgeqr2
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:132
dlarfg
subroutine dlarfg(N, ALPHA, X, INCX, TAU)
DLARFG generates an elementary reflector (Householder matrix).
Definition: dlarfg.f:108
xerbla
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
dlarf
subroutine dlarf(SIDE, M, N, V, INCV, TAU, C, LDC, WORK)
DLARF applies an elementary reflector to a general rectangular matrix.
Definition: dlarf.f:126