LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dgeqrt3.f
Go to the documentation of this file.
1*> \brief \b DGEQRT3 recursively computes a QR factorization of a general real or complex matrix using the compact WY representation of Q.
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download DGEQRT3 + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgeqrt3.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgeqrt3.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgeqrt3.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* RECURSIVE SUBROUTINE DGEQRT3( M, N, A, LDA, T, LDT, INFO )
20*
21* .. Scalar Arguments ..
22* INTEGER INFO, LDA, M, N, LDT
23* ..
24* .. Array Arguments ..
25* DOUBLE PRECISION A( LDA, * ), T( LDT, * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> DGEQRT3 recursively computes a QR factorization of a real M-by-N
35*> matrix A, using the compact WY representation of Q.
36*>
37*> Based on the algorithm of Elmroth and Gustavson,
38*> IBM J. Res. Develop. Vol 44 No. 4 July 2000.
39*> \endverbatim
40*
41* Arguments:
42* ==========
43*
44*> \param[in] M
45*> \verbatim
46*> M is INTEGER
47*> The number of rows of the matrix A. M >= N.
48*> \endverbatim
49*>
50*> \param[in] N
51*> \verbatim
52*> N is INTEGER
53*> The number of columns of the matrix A. N >= 0.
54*> \endverbatim
55*>
56*> \param[in,out] A
57*> \verbatim
58*> A is DOUBLE PRECISION array, dimension (LDA,N)
59*> On entry, the real M-by-N matrix A. On exit, the elements on and
60*> above the diagonal contain the N-by-N upper triangular matrix R; the
61*> elements below the diagonal are the columns of V. See below for
62*> further details.
63*> \endverbatim
64*>
65*> \param[in] LDA
66*> \verbatim
67*> LDA is INTEGER
68*> The leading dimension of the array A. LDA >= max(1,M).
69*> \endverbatim
70*>
71*> \param[out] T
72*> \verbatim
73*> T is DOUBLE PRECISION array, dimension (LDT,N)
74*> The N-by-N upper triangular factor of the block reflector.
75*> The elements on and above the diagonal contain the block
76*> reflector T; the elements below the diagonal are not used.
77*> See below for further details.
78*> \endverbatim
79*>
80*> \param[in] LDT
81*> \verbatim
82*> LDT is INTEGER
83*> The leading dimension of the array T. LDT >= max(1,N).
84*> \endverbatim
85*>
86*> \param[out] INFO
87*> \verbatim
88*> INFO is INTEGER
89*> = 0: successful exit
90*> < 0: if INFO = -i, the i-th argument had an illegal value
91*> \endverbatim
92*
93* Authors:
94* ========
95*
96*> \author Univ. of Tennessee
97*> \author Univ. of California Berkeley
98*> \author Univ. of Colorado Denver
99*> \author NAG Ltd.
100*
101*> \ingroup geqrt3
102*
103*> \par Further Details:
104* =====================
105*>
106*> \verbatim
107*>
108*> The matrix V stores the elementary reflectors H(i) in the i-th column
109*> below the diagonal. For example, if M=5 and N=3, the matrix V is
110*>
111*> V = ( 1 )
112*> ( v1 1 )
113*> ( v1 v2 1 )
114*> ( v1 v2 v3 )
115*> ( v1 v2 v3 )
116*>
117*> where the vi's represent the vectors which define H(i), which are returned
118*> in the matrix A. The 1's along the diagonal of V are not stored in A. The
119*> block reflector H is then given by
120*>
121*> H = I - V * T * V**T
122*>
123*> where V**T is the transpose of V.
124*>
125*> For details of the algorithm, see Elmroth and Gustavson (cited above).
126*> \endverbatim
127*>
128* =====================================================================
129 RECURSIVE SUBROUTINE dgeqrt3( M, N, A, LDA, T, LDT, INFO )
130*
131* -- LAPACK computational routine --
132* -- LAPACK is a software package provided by Univ. of Tennessee, --
133* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
134*
135* .. Scalar Arguments ..
136 INTEGER info, lda, m, n, ldt
137* ..
138* .. Array Arguments ..
139 DOUBLE PRECISION a( lda, * ), t( ldt, * )
140* ..
141*
142* =====================================================================
143*
144* .. Parameters ..
145 DOUBLE PRECISION one
146 parameter( one = 1.0d+00 )
147* ..
148* .. Local Scalars ..
149 INTEGER i, i1, j, j1, n1, n2, iinfo
150* ..
151* .. External Subroutines ..
152 EXTERNAL dlarfg, dtrmm, dgemm, xerbla
153* ..
154* .. Executable Statements ..
155*
156 info = 0
157 IF( n .LT. 0 ) THEN
158 info = -2
159 ELSE IF( m .LT. n ) THEN
160 info = -1
161 ELSE IF( lda .LT. max( 1, m ) ) THEN
162 info = -4
163 ELSE IF( ldt .LT. max( 1, n ) ) THEN
164 info = -6
165 END IF
166 IF( info.NE.0 ) THEN
167 CALL xerbla( 'DGEQRT3', -info )
168 RETURN
169 END IF
170*
171 IF( n.EQ.1 ) THEN
172*
173* Compute Householder transform when N=1
174*
175 CALL dlarfg( m, a(1,1), a( min( 2, m ), 1 ), 1, t(1,1) )
176*
177 ELSE
178*
179* Otherwise, split A into blocks...
180*
181 n1 = n/2
182 n2 = n-n1
183 j1 = min( n1+1, n )
184 i1 = min( n+1, m )
185*
186* Compute A(1:M,1:N1) <- (Y1,R1,T1), where Q1 = I - Y1 T1 Y1^H
187*
188 CALL dgeqrt3( m, n1, a, lda, t, ldt, iinfo )
189*
190* Compute A(1:M,J1:N) = Q1^H A(1:M,J1:N) [workspace: T(1:N1,J1:N)]
191*
192 DO j=1,n2
193 DO i=1,n1
194 t( i, j+n1 ) = a( i, j+n1 )
195 END DO
196 END DO
197 CALL dtrmm( 'L', 'L', 'T', 'U', n1, n2, one,
198 & a, lda, t( 1, j1 ), ldt )
199*
200 CALL dgemm( 'T', 'N', n1, n2, m-n1, one, a( j1, 1 ), lda,
201 & a( j1, j1 ), lda, one, t( 1, j1 ), ldt)
202*
203 CALL dtrmm( 'L', 'U', 'T', 'N', n1, n2, one,
204 & t, ldt, t( 1, j1 ), ldt )
205*
206 CALL dgemm( 'N', 'N', m-n1, n2, n1, -one, a( j1, 1 ), lda,
207 & t( 1, j1 ), ldt, one, a( j1, j1 ), lda )
208*
209 CALL dtrmm( 'L', 'L', 'N', 'U', n1, n2, one,
210 & a, lda, t( 1, j1 ), ldt )
211*
212 DO j=1,n2
213 DO i=1,n1
214 a( i, j+n1 ) = a( i, j+n1 ) - t( i, j+n1 )
215 END DO
216 END DO
217*
218* Compute A(J1:M,J1:N) <- (Y2,R2,T2) where Q2 = I - Y2 T2 Y2^H
219*
220 CALL dgeqrt3( m-n1, n2, a( j1, j1 ), lda,
221 & t( j1, j1 ), ldt, iinfo )
222*
223* Compute T3 = T(1:N1,J1:N) = -T1 Y1^H Y2 T2
224*
225 DO i=1,n1
226 DO j=1,n2
227 t( i, j+n1 ) = (a( j+n1, i ))
228 END DO
229 END DO
230*
231 CALL dtrmm( 'R', 'L', 'N', 'U', n1, n2, one,
232 & a( j1, j1 ), lda, t( 1, j1 ), ldt )
233*
234 CALL dgemm( 'T', 'N', n1, n2, m-n, one, a( i1, 1 ), lda,
235 & a( i1, j1 ), lda, one, t( 1, j1 ), ldt )
236*
237 CALL dtrmm( 'L', 'U', 'N', 'N', n1, n2, -one, t, ldt,
238 & t( 1, j1 ), ldt )
239*
240 CALL dtrmm( 'R', 'U', 'N', 'N', n1, n2, one,
241 & t( j1, j1 ), ldt, t( 1, j1 ), ldt )
242*
243* Y = (Y1,Y2); R = [ R1 A(1:N1,J1:N) ]; T = [T1 T3]
244* [ 0 R2 ] [ 0 T2]
245*
246 END IF
247*
248 RETURN
249*
250* End of DGEQRT3
251*
252 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
recursive subroutine dgeqrt3(m, n, a, lda, t, ldt, info)
DGEQRT3 recursively computes a QR factorization of a general real or complex matrix using the compact...
Definition dgeqrt3.f:130
subroutine dlarfg(n, alpha, x, incx, tau)
DLARFG generates an elementary reflector (Householder matrix).
Definition dlarfg.f:104
subroutine dtrmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
DTRMM
Definition dtrmm.f:177