LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
dgeqrfp.f
Go to the documentation of this file.
1 *> \brief \b DGEQRFP
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download DGEQRFP + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgeqrfp.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgeqrfp.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgeqrfp.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE DGEQRFP( M, N, A, LDA, TAU, WORK, LWORK, INFO )
22 *
23 * .. Scalar Arguments ..
24 * INTEGER INFO, LDA, LWORK, 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 *> DGEQRFP computes a QR factorization of a real M-by-N matrix A:
37 *> A = Q * R. The diagonal entries of R are nonnegative.
38 *> \endverbatim
39 *
40 * Arguments:
41 * ==========
42 *
43 *> \param[in] M
44 *> \verbatim
45 *> M is INTEGER
46 *> The number of rows of the matrix A. M >= 0.
47 *> \endverbatim
48 *>
49 *> \param[in] N
50 *> \verbatim
51 *> N is INTEGER
52 *> The number of columns of the matrix A. N >= 0.
53 *> \endverbatim
54 *>
55 *> \param[in,out] A
56 *> \verbatim
57 *> A is DOUBLE PRECISION array, dimension (LDA,N)
58 *> On entry, the M-by-N matrix A.
59 *> On exit, the elements on and above the diagonal of the array
60 *> contain the min(M,N)-by-N upper trapezoidal matrix R (R is
61 *> upper triangular if m >= n). The diagonal entries of R
62 *> are nonnegative; the elements below the diagonal,
63 *> with the array TAU, represent the orthogonal matrix Q as a
64 *> product of min(m,n) elementary reflectors (see Further
65 *> Details).
66 *> \endverbatim
67 *>
68 *> \param[in] LDA
69 *> \verbatim
70 *> LDA is INTEGER
71 *> The leading dimension of the array A. LDA >= max(1,M).
72 *> \endverbatim
73 *>
74 *> \param[out] TAU
75 *> \verbatim
76 *> TAU is DOUBLE PRECISION array, dimension (min(M,N))
77 *> The scalar factors of the elementary reflectors (see Further
78 *> Details).
79 *> \endverbatim
80 *>
81 *> \param[out] WORK
82 *> \verbatim
83 *> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
84 *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
85 *> \endverbatim
86 *>
87 *> \param[in] LWORK
88 *> \verbatim
89 *> LWORK is INTEGER
90 *> The dimension of the array WORK. LWORK >= max(1,N).
91 *> For optimum performance LWORK >= N*NB, where NB is
92 *> the optimal blocksize.
93 *>
94 *> If LWORK = -1, then a workspace query is assumed; the routine
95 *> only calculates the optimal size of the WORK array, returns
96 *> this value as the first entry of the WORK array, and no error
97 *> message related to LWORK is issued by XERBLA.
98 *> \endverbatim
99 *>
100 *> \param[out] INFO
101 *> \verbatim
102 *> INFO is INTEGER
103 *> = 0: successful exit
104 *> < 0: if INFO = -i, the i-th argument had an illegal value
105 *> \endverbatim
106 *
107 * Authors:
108 * ========
109 *
110 *> \author Univ. of Tennessee
111 *> \author Univ. of California Berkeley
112 *> \author Univ. of Colorado Denver
113 *> \author NAG Ltd.
114 *
115 *> \date November 2015
116 *
117 *> \ingroup doubleGEcomputational
118 *
119 *> \par Further Details:
120 * =====================
121 *>
122 *> \verbatim
123 *>
124 *> The matrix Q is represented as a product of elementary reflectors
125 *>
126 *> Q = H(1) H(2) . . . H(k), where k = min(m,n).
127 *>
128 *> Each H(i) has the form
129 *>
130 *> H(i) = I - tau * v * v**T
131 *>
132 *> where tau is a real scalar, and v is a real vector with
133 *> v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i),
134 *> and tau in TAU(i).
135 *>
136 *> See Lapack Working Note 203 for details
137 *> \endverbatim
138 *>
139 * =====================================================================
140  SUBROUTINE dgeqrfp( M, N, A, LDA, TAU, WORK, LWORK, INFO )
141 *
142 * -- LAPACK computational routine (version 3.6.0) --
143 * -- LAPACK is a software package provided by Univ. of Tennessee, --
144 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
145 * November 2015
146 *
147 * .. Scalar Arguments ..
148  INTEGER INFO, LDA, LWORK, M, N
149 * ..
150 * .. Array Arguments ..
151  DOUBLE PRECISION A( lda, * ), TAU( * ), WORK( * )
152 * ..
153 *
154 * =====================================================================
155 *
156 * .. Local Scalars ..
157  LOGICAL LQUERY
158  INTEGER I, IB, IINFO, IWS, K, LDWORK, LWKOPT, NB,
159  $ nbmin, nx
160 * ..
161 * .. External Subroutines ..
162  EXTERNAL dgeqr2p, dlarfb, dlarft, xerbla
163 * ..
164 * .. Intrinsic Functions ..
165  INTRINSIC max, min
166 * ..
167 * .. External Functions ..
168  INTEGER ILAENV
169  EXTERNAL ilaenv
170 * ..
171 * .. Executable Statements ..
172 *
173 * Test the input arguments
174 *
175  info = 0
176  nb = ilaenv( 1, 'DGEQRF', ' ', m, n, -1, -1 )
177  lwkopt = n*nb
178  work( 1 ) = lwkopt
179  lquery = ( lwork.EQ.-1 )
180  IF( m.LT.0 ) THEN
181  info = -1
182  ELSE IF( n.LT.0 ) THEN
183  info = -2
184  ELSE IF( lda.LT.max( 1, m ) ) THEN
185  info = -4
186  ELSE IF( lwork.LT.max( 1, n ) .AND. .NOT.lquery ) THEN
187  info = -7
188  END IF
189  IF( info.NE.0 ) THEN
190  CALL xerbla( 'DGEQRFP', -info )
191  RETURN
192  ELSE IF( lquery ) THEN
193  RETURN
194  END IF
195 *
196 * Quick return if possible
197 *
198  k = min( m, n )
199  IF( k.EQ.0 ) THEN
200  work( 1 ) = 1
201  RETURN
202  END IF
203 *
204  nbmin = 2
205  nx = 0
206  iws = n
207  IF( nb.GT.1 .AND. nb.LT.k ) THEN
208 *
209 * Determine when to cross over from blocked to unblocked code.
210 *
211  nx = max( 0, ilaenv( 3, 'DGEQRF', ' ', m, n, -1, -1 ) )
212  IF( nx.LT.k ) THEN
213 *
214 * Determine if workspace is large enough for blocked code.
215 *
216  ldwork = n
217  iws = ldwork*nb
218  IF( lwork.LT.iws ) THEN
219 *
220 * Not enough workspace to use optimal NB: reduce NB and
221 * determine the minimum value of NB.
222 *
223  nb = lwork / ldwork
224  nbmin = max( 2, ilaenv( 2, 'DGEQRF', ' ', m, n, -1,
225  $ -1 ) )
226  END IF
227  END IF
228  END IF
229 *
230  IF( nb.GE.nbmin .AND. nb.LT.k .AND. nx.LT.k ) THEN
231 *
232 * Use blocked code initially
233 *
234  DO 10 i = 1, k - nx, nb
235  ib = min( k-i+1, nb )
236 *
237 * Compute the QR factorization of the current block
238 * A(i:m,i:i+ib-1)
239 *
240  CALL dgeqr2p( m-i+1, ib, a( i, i ), lda, tau( i ), work,
241  $ iinfo )
242  IF( i+ib.LE.n ) THEN
243 *
244 * Form the triangular factor of the block reflector
245 * H = H(i) H(i+1) . . . H(i+ib-1)
246 *
247  CALL dlarft( 'Forward', 'Columnwise', m-i+1, ib,
248  $ a( i, i ), lda, tau( i ), work, ldwork )
249 *
250 * Apply H**T to A(i:m,i+ib:n) from the left
251 *
252  CALL dlarfb( 'Left', 'Transpose', 'Forward',
253  $ 'Columnwise', m-i+1, n-i-ib+1, ib,
254  $ a( i, i ), lda, work, ldwork, a( i, i+ib ),
255  $ lda, work( ib+1 ), ldwork )
256  END IF
257  10 CONTINUE
258  ELSE
259  i = 1
260  END IF
261 *
262 * Use unblocked code to factor the last or only block.
263 *
264  IF( i.LE.k )
265  $ CALL dgeqr2p( m-i+1, n-i+1, a( i, i ), lda, tau( i ), work,
266  $ iinfo )
267 *
268  work( 1 ) = iws
269  RETURN
270 *
271 * End of DGEQRFP
272 *
273  END
subroutine dlarfb(SIDE, TRANS, DIRECT, STOREV, M, N, K, V, LDV, T, LDT, C, LDC, WORK, LDWORK)
DLARFB applies a block reflector or its transpose to a general rectangular matrix.
Definition: dlarfb.f:197
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dlarft(DIRECT, STOREV, N, K, V, LDV, TAU, T, LDT)
DLARFT forms the triangular factor T of a block reflector H = I - vtvH
Definition: dlarft.f:165
subroutine dgeqr2p(M, N, A, LDA, TAU, WORK, INFO)
DGEQR2P computes the QR factorization of a general rectangular matrix with non-negative diagonal elem...
Definition: dgeqr2p.f:126
subroutine dgeqrfp(M, N, A, LDA, TAU, WORK, LWORK, INFO)
DGEQRFP
Definition: dgeqrfp.f:141