LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zgeqp3rk.f
Go to the documentation of this file.
1*> \brief \b ZGEQP3RK computes a truncated Householder QR factorization with column pivoting of a complex m-by-n matrix A by using Level 3 BLAS and overwrites m-by-nrhs matrix B with Q**H * B.
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download ZGEQP3RK + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zgeqp3rk.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zgeqp3rk.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zgeqp3rk.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE ZGEQP3RK( M, N, NRHS, KMAX, ABSTOL, RELTOL, A, LDA,
20* $ K, MAXC2NRMK, RELMAXC2NRMK, JPIV, TAU,
21* $ WORK, LWORK, RWORK, IWORK, INFO )
22* IMPLICIT NONE
23*
24* .. Scalar Arguments ..
25* INTEGER INFO, K, KMAX, LDA, LWORK, M, N, NRHS
26* DOUBLE PRECISION ABSTOL, MAXC2NRMK, RELMAXC2NRMK, RELTOL
27* ..
28* .. Array Arguments ..
29* INTEGER IWORK( * ), JPIV( * )
30* DOUBLE PRECISION RWORK( * )
31* COMPLEX*16 A( LDA, * ), TAU( * ), WORK( * )
32* ..
33*
34*
35*> \par Purpose:
36* =============
37*>
38*> \verbatim
39*>
40*> ZGEQP3RK performs two tasks simultaneously:
41*>
42*> Task 1: The routine computes a truncated (rank K) or full rank
43*> Householder QR factorization with column pivoting of a complex
44*> M-by-N matrix A using Level 3 BLAS. K is the number of columns
45*> that were factorized, i.e. factorization rank of the
46*> factor R, K <= min(M,N).
47*>
48*> A * P(K) = Q(K) * R(K) =
49*>
50*> = Q(K) * ( R11(K) R12(K) ) = Q(K) * ( R(K)_approx )
51*> ( 0 R22(K) ) ( 0 R(K)_residual ),
52*>
53*> where:
54*>
55*> P(K) is an N-by-N permutation matrix;
56*> Q(K) is an M-by-M unitary matrix;
57*> R(K)_approx = ( R11(K), R12(K) ) is a rank K approximation of the
58*> full rank factor R with K-by-K upper-triangular
59*> R11(K) and K-by-N rectangular R12(K). The diagonal
60*> entries of R11(K) appear in non-increasing order
61*> of absolute value, and absolute values of all of
62*> them exceed the maximum column 2-norm of R22(K)
63*> up to roundoff error.
64*> R(K)_residual = R22(K) is the residual of a rank K approximation
65*> of the full rank factor R. It is a
66*> an (M-K)-by-(N-K) rectangular matrix;
67*> 0 is a an (M-K)-by-K zero matrix.
68*>
69*> Task 2: At the same time, the routine overwrites a complex M-by-NRHS
70*> matrix B with Q(K)**H * B using Level 3 BLAS.
71*>
72*> =====================================================================
73*>
74*> The matrices A and B are stored on input in the array A as
75*> the left and right blocks A(1:M,1:N) and A(1:M, N+1:N+NRHS)
76*> respectively.
77*>
78*> N NRHS
79*> array_A = M [ mat_A, mat_B ]
80*>
81*> The truncation criteria (i.e. when to stop the factorization)
82*> can be any of the following:
83*>
84*> 1) The input parameter KMAX, the maximum number of columns
85*> KMAX to factorize, i.e. the factorization rank is limited
86*> to KMAX. If KMAX >= min(M,N), the criterion is not used.
87*>
88*> 2) The input parameter ABSTOL, the absolute tolerance for
89*> the maximum column 2-norm of the residual matrix R22(K). This
90*> means that the factorization stops if this norm is less or
91*> equal to ABSTOL. If ABSTOL < 0.0, the criterion is not used.
92*>
93*> 3) The input parameter RELTOL, the tolerance for the maximum
94*> column 2-norm matrix of the residual matrix R22(K) divided
95*> by the maximum column 2-norm of the original matrix A, which
96*> is equal to abs(R(1,1)). This means that the factorization stops
97*> when the ratio of the maximum column 2-norm of R22(K) to
98*> the maximum column 2-norm of A is less than or equal to RELTOL.
99*> If RELTOL < 0.0, the criterion is not used.
100*>
101*> 4) In case both stopping criteria ABSTOL or RELTOL are not used,
102*> and when the residual matrix R22(K) is a zero matrix in some
103*> factorization step K. ( This stopping criterion is implicit. )
104*>
105*> The algorithm stops when any of these conditions is first
106*> satisfied, otherwise the whole matrix A is factorized.
107*>
108*> To factorize the whole matrix A, use the values
109*> KMAX >= min(M,N), ABSTOL < 0.0 and RELTOL < 0.0.
110*>
111*> The routine returns:
112*> a) Q(K), R(K)_approx = ( R11(K), R12(K) ),
113*> R(K)_residual = R22(K), P(K), i.e. the resulting matrices
114*> of the factorization; P(K) is represented by JPIV,
115*> ( if K = min(M,N), R(K)_approx is the full factor R,
116*> and there is no residual matrix R(K)_residual);
117*> b) K, the number of columns that were factorized,
118*> i.e. factorization rank;
119*> c) MAXC2NRMK, the maximum column 2-norm of the residual
120*> matrix R(K)_residual = R22(K),
121*> ( if K = min(M,N), MAXC2NRMK = 0.0 );
122*> d) RELMAXC2NRMK equals MAXC2NRMK divided by MAXC2NRM, the maximum
123*> column 2-norm of the original matrix A, which is equal
124*> to abs(R(1,1)), ( if K = min(M,N), RELMAXC2NRMK = 0.0 );
125*> e) Q(K)**H * B, the matrix B with the unitary
126*> transformation Q(K)**H applied on the left.
127*>
128*> The N-by-N permutation matrix P(K) is stored in a compact form in
129*> the integer array JPIV. For 1 <= j <= N, column j
130*> of the matrix A was interchanged with column JPIV(j).
131*>
132*> The M-by-M unitary matrix Q is represented as a product
133*> of elementary Householder reflectors
134*>
135*> Q(K) = H(1) * H(2) * . . . * H(K),
136*>
137*> where K is the number of columns that were factorized.
138*>
139*> Each H(j) has the form
140*>
141*> H(j) = I - tau * v * v**H,
142*>
143*> where 1 <= j <= K and
144*> I is an M-by-M identity matrix,
145*> tau is a complex scalar,
146*> v is a complex vector with v(1:j-1) = 0 and v(j) = 1.
147*>
148*> v(j+1:M) is stored on exit in A(j+1:M,j) and tau in TAU(j).
149*>
150*> See the Further Details section for more information.
151*> \endverbatim
152*
153* Arguments:
154* ==========
155*
156*> \param[in] M
157*> \verbatim
158*> M is INTEGER
159*> The number of rows of the matrix A. M >= 0.
160*> \endverbatim
161*>
162*> \param[in] N
163*> \verbatim
164*> N is INTEGER
165*> The number of columns of the matrix A. N >= 0.
166*> \endverbatim
167*>
168*> \param[in] NRHS
169*> \verbatim
170*> NRHS is INTEGER
171*> The number of right hand sides, i.e. the number of
172*> columns of the matrix B. NRHS >= 0.
173*> \endverbatim
174*>
175*> \param[in] KMAX
176*> \verbatim
177*> KMAX is INTEGER
178*>
179*> The first factorization stopping criterion. KMAX >= 0.
180*>
181*> The maximum number of columns of the matrix A to factorize,
182*> i.e. the maximum factorization rank.
183*>
184*> a) If KMAX >= min(M,N), then this stopping criterion
185*> is not used, the routine factorizes columns
186*> depending on ABSTOL and RELTOL.
187*>
188*> b) If KMAX = 0, then this stopping criterion is
189*> satisfied on input and the routine exits immediately.
190*> This means that the factorization is not performed,
191*> the matrices A and B are not modified, and
192*> the matrix A is itself the residual.
193*> \endverbatim
194*>
195*> \param[in] ABSTOL
196*> \verbatim
197*> ABSTOL is DOUBLE PRECISION
198*>
199*> The second factorization stopping criterion, cannot be NaN.
200*>
201*> The absolute tolerance (stopping threshold) for
202*> maximum column 2-norm of the residual matrix R22(K).
203*> The algorithm converges (stops the factorization) when
204*> the maximum column 2-norm of the residual matrix R22(K)
205*> is less than or equal to ABSTOL. Let SAFMIN = DLAMCH('S').
206*>
207*> a) If ABSTOL is NaN, then no computation is performed
208*> and an error message ( INFO = -5 ) is issued
209*> by XERBLA.
210*>
211*> b) If ABSTOL < 0.0, then this stopping criterion is not
212*> used, the routine factorizes columns depending
213*> on KMAX and RELTOL.
214*> This includes the case ABSTOL = -Inf.
215*>
216*> c) If 0.0 <= ABSTOL < 2*SAFMIN, then ABSTOL = 2*SAFMIN
217*> is used. This includes the case ABSTOL = -0.0.
218*>
219*> d) If 2*SAFMIN <= ABSTOL then the input value
220*> of ABSTOL is used.
221*>
222*> Let MAXC2NRM be the maximum column 2-norm of the
223*> whole original matrix A.
224*> If ABSTOL chosen above is >= MAXC2NRM, then this
225*> stopping criterion is satisfied on input and routine exits
226*> immediately after MAXC2NRM is computed. The routine
227*> returns MAXC2NRM in MAXC2NORMK,
228*> and 1.0 in RELMAXC2NORMK.
229*> This includes the case ABSTOL = +Inf. This means that the
230*> factorization is not performed, the matrices A and B are not
231*> modified, and the matrix A is itself the residual.
232*> \endverbatim
233*>
234*> \param[in] RELTOL
235*> \verbatim
236*> RELTOL is DOUBLE PRECISION
237*>
238*> The third factorization stopping criterion, cannot be NaN.
239*>
240*> The tolerance (stopping threshold) for the ratio
241*> abs(R(K+1,K+1))/abs(R(1,1)) of the maximum column 2-norm of
242*> the residual matrix R22(K) to the maximum column 2-norm of
243*> the original matrix A. The algorithm converges (stops the
244*> factorization), when abs(R(K+1,K+1))/abs(R(1,1)) A is less
245*> than or equal to RELTOL. Let EPS = DLAMCH('E').
246*>
247*> a) If RELTOL is NaN, then no computation is performed
248*> and an error message ( INFO = -6 ) is issued
249*> by XERBLA.
250*>
251*> b) If RELTOL < 0.0, then this stopping criterion is not
252*> used, the routine factorizes columns depending
253*> on KMAX and ABSTOL.
254*> This includes the case RELTOL = -Inf.
255*>
256*> c) If 0.0 <= RELTOL < EPS, then RELTOL = EPS is used.
257*> This includes the case RELTOL = -0.0.
258*>
259*> d) If EPS <= RELTOL then the input value of RELTOL
260*> is used.
261*>
262*> Let MAXC2NRM be the maximum column 2-norm of the
263*> whole original matrix A.
264*> If RELTOL chosen above is >= 1.0, then this stopping
265*> criterion is satisfied on input and routine exits
266*> immediately after MAXC2NRM is computed.
267*> The routine returns MAXC2NRM in MAXC2NORMK,
268*> and 1.0 in RELMAXC2NORMK.
269*> This includes the case RELTOL = +Inf. This means that the
270*> factorization is not performed, the matrices A and B are not
271*> modified, and the matrix A is itself the residual.
272*>
273*> NOTE: We recommend that RELTOL satisfy
274*> min( 10*max(M,N)*EPS, sqrt(EPS) ) <= RELTOL
275*> \endverbatim
276*>
277*> \param[in,out] A
278*> \verbatim
279*> A is COMPLEX*16 array, dimension (LDA,N+NRHS)
280*>
281*> On entry:
282*>
283*> a) The subarray A(1:M,1:N) contains the M-by-N matrix A.
284*> b) The subarray A(1:M,N+1:N+NRHS) contains the M-by-NRHS
285*> matrix B.
286*>
287*> N NRHS
288*> array_A = M [ mat_A, mat_B ]
289*>
290*> On exit:
291*>
292*> a) The subarray A(1:M,1:N) contains parts of the factors
293*> of the matrix A:
294*>
295*> 1) If K = 0, A(1:M,1:N) contains the original matrix A.
296*> 2) If K > 0, A(1:M,1:N) contains parts of the
297*> factors:
298*>
299*> 1. The elements below the diagonal of the subarray
300*> A(1:M,1:K) together with TAU(1:K) represent the
301*> unitary matrix Q(K) as a product of K Householder
302*> elementary reflectors.
303*>
304*> 2. The elements on and above the diagonal of
305*> the subarray A(1:K,1:N) contain K-by-N
306*> upper-trapezoidal matrix
307*> R(K)_approx = ( R11(K), R12(K) ).
308*> NOTE: If K=min(M,N), i.e. full rank factorization,
309*> then R_approx(K) is the full factor R which
310*> is upper-trapezoidal. If, in addition, M>=N,
311*> then R is upper-triangular.
312*>
313*> 3. The subarray A(K+1:M,K+1:N) contains (M-K)-by-(N-K)
314*> rectangular matrix R(K)_residual = R22(K).
315*>
316*> b) If NRHS > 0, the subarray A(1:M,N+1:N+NRHS) contains
317*> the M-by-NRHS product Q(K)**H * B.
318*> \endverbatim
319*>
320*> \param[in] LDA
321*> \verbatim
322*> LDA is INTEGER
323*> The leading dimension of the array A. LDA >= max(1,M).
324*> This is the leading dimension for both matrices, A and B.
325*> \endverbatim
326*>
327*> \param[out] K
328*> \verbatim
329*> K is INTEGER
330*> Factorization rank of the matrix A, i.e. the rank of
331*> the factor R, which is the same as the number of non-zero
332*> rows of the factor R. 0 <= K <= min(M,KMAX,N).
333*>
334*> K also represents the number of non-zero Householder
335*> vectors.
336*>
337*> NOTE: If K = 0, a) the arrays A and B are not modified;
338*> b) the array TAU(1:min(M,N)) is set to ZERO,
339*> if the matrix A does not contain NaN,
340*> otherwise the elements TAU(1:min(M,N))
341*> are undefined;
342*> c) the elements of the array JPIV are set
343*> as follows: for j = 1:N, JPIV(j) = j.
344*> \endverbatim
345*>
346*> \param[out] MAXC2NRMK
347*> \verbatim
348*> MAXC2NRMK is DOUBLE PRECISION
349*> The maximum column 2-norm of the residual matrix R22(K),
350*> when the factorization stopped at rank K. MAXC2NRMK >= 0.
351*>
352*> a) If K = 0, i.e. the factorization was not performed,
353*> the matrix A was not modified and is itself a residual
354*> matrix, then MAXC2NRMK equals the maximum column 2-norm
355*> of the original matrix A.
356*>
357*> b) If 0 < K < min(M,N), then MAXC2NRMK is returned.
358*>
359*> c) If K = min(M,N), i.e. the whole matrix A was
360*> factorized and there is no residual matrix,
361*> then MAXC2NRMK = 0.0.
362*>
363*> NOTE: MAXC2NRMK in the factorization step K would equal
364*> R(K+1,K+1) in the next factorization step K+1.
365*> \endverbatim
366*>
367*> \param[out] RELMAXC2NRMK
368*> \verbatim
369*> RELMAXC2NRMK is DOUBLE PRECISION
370*> The ratio MAXC2NRMK / MAXC2NRM of the maximum column
371*> 2-norm of the residual matrix R22(K) (when the factorization
372*> stopped at rank K) to the maximum column 2-norm of the
373*> whole original matrix A. RELMAXC2NRMK >= 0.
374*>
375*> a) If K = 0, i.e. the factorization was not performed,
376*> the matrix A was not modified and is itself a residual
377*> matrix, then RELMAXC2NRMK = 1.0.
378*>
379*> b) If 0 < K < min(M,N), then
380*> RELMAXC2NRMK = MAXC2NRMK / MAXC2NRM is returned.
381*>
382*> c) If K = min(M,N), i.e. the whole matrix A was
383*> factorized and there is no residual matrix,
384*> then RELMAXC2NRMK = 0.0.
385*>
386*> NOTE: RELMAXC2NRMK in the factorization step K would equal
387*> abs(R(K+1,K+1))/abs(R(1,1)) in the next factorization
388*> step K+1.
389*> \endverbatim
390*>
391*> \param[out] JPIV
392*> \verbatim
393*> JPIV is INTEGER array, dimension (N)
394*> Column pivot indices. For 1 <= j <= N, column j
395*> of the matrix A was interchanged with column JPIV(j).
396*>
397*> The elements of the array JPIV(1:N) are always set
398*> by the routine, for example, even when no columns
399*> were factorized, i.e. when K = 0, the elements are
400*> set as JPIV(j) = j for j = 1:N.
401*> \endverbatim
402*>
403*> \param[out] TAU
404*> \verbatim
405*> TAU is COMPLEX*16 array, dimension (min(M,N))
406*> The scalar factors of the elementary reflectors.
407*>
408*> If 0 < K <= min(M,N), only the elements TAU(1:K) of
409*> the array TAU are modified by the factorization.
410*> After the factorization computed, if no NaN was found
411*> during the factorization, the remaining elements
412*> TAU(K+1:min(M,N)) are set to zero, otherwise the
413*> elements TAU(K+1:min(M,N)) are not set and therefore
414*> undefined.
415*> ( If K = 0, all elements of TAU are set to zero, if
416*> the matrix A does not contain NaN. )
417*> \endverbatim
418*>
419*> \param[out] WORK
420*> \verbatim
421*> WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
422*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
423*> \endverbatim
424*>
425*> \param[in] LWORK
426*> \verbatim
427*> LWORK is INTEGER
428*> The dimension of the array WORK.
429*> LWORK >= 1, if MIN(M,N) = 0, and
430*> LWORK >= N+NRHS-1, otherwise.
431*> For optimal performance LWORK >= NB*( N+NRHS+1 ),
432*> where NB is the optimal block size for ZGEQP3RK returned
433*> by ILAENV. Minimal block size MINNB=2.
434*>
435*> NOTE: The decision, whether to use unblocked BLAS 2
436*> or blocked BLAS 3 code is based not only on the dimension
437*> LWORK of the availbale workspace WORK, but also also on the
438*> matrix A dimension N via crossover point NX returned
439*> by ILAENV. (For N less than NX, unblocked code should be
440*> used.)
441*>
442*> If LWORK = -1, then a workspace query is assumed;
443*> the routine only calculates the optimal size of the WORK
444*> array, returns this value as the first entry of the WORK
445*> array, and no error message related to LWORK is issued
446*> by XERBLA.
447*> \endverbatim
448*>
449*> \param[out] RWORK
450*> \verbatim
451*> RWORK is DOUBLE PRECISION array, dimension (2*N)
452*> \endverbatim
453*>
454*> \param[out] IWORK
455*> \verbatim
456*> IWORK is INTEGER array, dimension (N-1).
457*> Is a work array. ( IWORK is used to store indices
458*> of "bad" columns for norm downdating in the residual
459*> matrix in the blocked step auxiliary subroutine ZLAQP3RK ).
460*> \endverbatim
461*>
462*> \param[out] INFO
463*> \verbatim
464*> INFO is INTEGER
465*> 1) INFO = 0: successful exit.
466*> 2) INFO < 0: if INFO = -i, the i-th argument had an
467*> illegal value.
468*> 3) If INFO = j_1, where 1 <= j_1 <= N, then NaN was
469*> detected and the routine stops the computation.
470*> The j_1-th column of the matrix A or the j_1-th
471*> element of array TAU contains the first occurrence
472*> of NaN in the factorization step K+1 ( when K columns
473*> have been factorized ).
474*>
475*> On exit:
476*> K is set to the number of
477*> factorized columns without
478*> exception.
479*> MAXC2NRMK is set to NaN.
480*> RELMAXC2NRMK is set to NaN.
481*> TAU(K+1:min(M,N)) is not set and contains undefined
482*> elements. If j_1=K+1, TAU(K+1)
483*> may contain NaN.
484*> 4) If INFO = j_2, where N+1 <= j_2 <= 2*N, then no NaN
485*> was detected, but +Inf (or -Inf) was detected and
486*> the routine continues the computation until completion.
487*> The (j_2-N)-th column of the matrix A contains the first
488*> occurrence of +Inf (or -Inf) in the factorization
489*> step K+1 ( when K columns have been factorized ).
490*> \endverbatim
491*
492* Authors:
493* ========
494*
495*> \author Univ. of Tennessee
496*> \author Univ. of California Berkeley
497*> \author Univ. of Colorado Denver
498*> \author NAG Ltd.
499*
500*> \ingroup geqp3rk
501*
502*> \par Further Details:
503* =====================
504*
505*> \verbatim
506*> ZGEQP3RK is based on the same BLAS3 Householder QR factorization
507*> algorithm with column pivoting as in ZGEQP3 routine which uses
508*> ZLARFG routine to generate Householder reflectors
509*> for QR factorization.
510*>
511*> We can also write:
512*>
513*> A = A_approx(K) + A_residual(K)
514*>
515*> The low rank approximation matrix A(K)_approx from
516*> the truncated QR factorization of rank K of the matrix A is:
517*>
518*> A(K)_approx = Q(K) * ( R(K)_approx ) * P(K)**T
519*> ( 0 0 )
520*>
521*> = Q(K) * ( R11(K) R12(K) ) * P(K)**T
522*> ( 0 0 )
523*>
524*> The residual A_residual(K) of the matrix A is:
525*>
526*> A_residual(K) = Q(K) * ( 0 0 ) * P(K)**T =
527*> ( 0 R(K)_residual )
528*>
529*> = Q(K) * ( 0 0 ) * P(K)**T
530*> ( 0 R22(K) )
531*>
532*> The truncated (rank K) factorization guarantees that
533*> the maximum column 2-norm of A_residual(K) is less than
534*> or equal to MAXC2NRMK up to roundoff error.
535*>
536*> NOTE: An approximation of the null vectors
537*> of A can be easily computed from R11(K)
538*> and R12(K):
539*>
540*> Null( A(K) )_approx = P * ( inv(R11(K)) * R12(K) )
541*> ( -I )
542*>
543*> \endverbatim
544*
545*> \par References:
546* ================
547*> [1] A Level 3 BLAS QR factorization algorithm with column pivoting developed in 1996.
548*> G. Quintana-Orti, Depto. de Informatica, Universidad Jaime I, Spain.
549*> X. Sun, Computer Science Dept., Duke University, USA.
550*> C. H. Bischof, Math. and Comp. Sci. Div., Argonne National Lab, USA.
551*> A BLAS-3 version of the QR factorization with column pivoting.
552*> LAPACK Working Note 114
553*> <a href="https://www.netlib.org/lapack/lawnspdf/lawn114.pdf">https://www.netlib.org/lapack/lawnspdf/lawn114.pdf</a>
554*> and in
555*> SIAM J. Sci. Comput., 19(5):1486-1494, Sept. 1998.
556*> <a href="https://doi.org/10.1137/S1064827595296732">https://doi.org/10.1137/S1064827595296732</a>
557*>
558*> [2] A partial column norm updating strategy developed in 2006.
559*> Z. Drmac and Z. Bujanovic, Dept. of Math., University of Zagreb, Croatia.
560*> On the failure of rank revealing QR factorization software – a case study.
561*> LAPACK Working Note 176.
562*> <a href="http://www.netlib.org/lapack/lawnspdf/lawn176.pdf">http://www.netlib.org/lapack/lawnspdf/lawn176.pdf</a>
563*> and in
564*> ACM Trans. Math. Softw. 35, 2, Article 12 (July 2008), 28 pages.
565*> <a href="https://doi.org/10.1145/1377612.1377616">https://doi.org/10.1145/1377612.1377616</a>
566*
567*> \par Contributors:
568* ==================
569*>
570*> \verbatim
571*>
572*> November 2023, Igor Kozachenko, James Demmel,
573*> EECS Department,
574*> University of California, Berkeley, USA.
575*>
576*> \endverbatim
577*
578* =====================================================================
579 SUBROUTINE zgeqp3rk( M, N, NRHS, KMAX, ABSTOL, RELTOL, A, LDA,
580 $ K, MAXC2NRMK, RELMAXC2NRMK, JPIV, TAU,
581 $ WORK, LWORK, RWORK, IWORK, INFO )
582 IMPLICIT NONE
583*
584* -- LAPACK computational routine --
585* -- LAPACK is a software package provided by Univ. of Tennessee, --
586* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
587*
588* .. Scalar Arguments ..
589 INTEGER INFO, K, KF, KMAX, LDA, LWORK, M, N, NRHS
590 DOUBLE PRECISION ABSTOL, MAXC2NRMK, RELMAXC2NRMK, RELTOL
591* ..
592* .. Array Arguments ..
593 INTEGER IWORK( * ), JPIV( * )
594 DOUBLE PRECISION RWORK( * )
595 COMPLEX*16 A( LDA, * ), TAU( * ), WORK( * )
596* ..
597*
598* =====================================================================
599*
600* .. Parameters ..
601 INTEGER INB, INBMIN, IXOVER
602 PARAMETER ( INB = 1, inbmin = 2, ixover = 3 )
603 DOUBLE PRECISION ZERO, ONE, TWO
604 parameter( zero = 0.0d+0, one = 1.0d+0, two = 2.0d+0 )
605 COMPLEX*16 CZERO
606 parameter( czero = ( 0.0d+0, 0.0d+0 ) )
607* ..
608* .. Local Scalars ..
609 LOGICAL LQUERY, DONE
610 INTEGER IINFO, IOFFSET, IWS, J, JB, JBF, JMAXB, JMAX,
611 $ jmaxc2nrm, kp1, lwkopt, minmn, n_sub, nb,
612 $ nbmin, nx
613 DOUBLE PRECISION EPS, HUGEVAL, MAXC2NRM, SAFMIN
614* ..
615* .. External Subroutines ..
616 EXTERNAL zlaqp2rk, zlaqp3rk, xerbla
617* ..
618* .. External Functions ..
619 LOGICAL DISNAN
620 INTEGER IDAMAX, ILAENV
621 DOUBLE PRECISION DLAMCH, DZNRM2
622 EXTERNAL disnan, dlamch, dznrm2, idamax, ilaenv
623* ..
624* .. Intrinsic Functions ..
625 INTRINSIC dcmplx, max, min
626* ..
627* .. Executable Statements ..
628*
629* Test input arguments
630* ====================
631*
632 info = 0
633 lquery = ( lwork.EQ.-1 )
634 IF( m.LT.0 ) THEN
635 info = -1
636 ELSE IF( n.LT.0 ) THEN
637 info = -2
638 ELSE IF( nrhs.LT.0 ) THEN
639 info = -3
640 ELSE IF( kmax.LT.0 ) THEN
641 info = -4
642 ELSE IF( disnan( abstol ) ) THEN
643 info = -5
644 ELSE IF( disnan( reltol ) ) THEN
645 info = -6
646 ELSE IF( lda.LT.max( 1, m ) ) THEN
647 info = -8
648 END IF
649*
650* If the input parameters M, N, NRHS, KMAX, LDA are valid:
651* a) Test the input workspace size LWORK for the minimum
652* size requirement IWS.
653* b) Determine the optimal block size NB and optimal
654* workspace size LWKOPT to be returned in WORK(1)
655* in case of (1) LWORK < IWS, (2) LQUERY = .TRUE.,
656* (3) when routine exits.
657* Here, IWS is the miminum workspace required for unblocked
658* code.
659*
660 IF( info.EQ.0 ) THEN
661 minmn = min( m, n )
662 IF( minmn.EQ.0 ) THEN
663 iws = 1
664 lwkopt = 1
665 ELSE
666*
667* Minimal workspace size in case of using only unblocked
668* BLAS 2 code in ZLAQP2RK.
669* 1) ZLAQP2RK: N+NRHS-1 to use in WORK array that is used
670* in ZLARF1F subroutine inside ZLAQP2RK to apply an
671* elementary reflector from the left.
672* TOTAL_WORK_SIZE = 3*N + NRHS - 1
673*
674 iws = n + nrhs - 1
675*
676* Assign to NB optimal block size.
677*
678 nb = ilaenv( inb, 'ZGEQP3RK', ' ', m, n, -1, -1 )
679*
680* A formula for the optimal workspace size in case of using
681* both unblocked BLAS 2 in ZLAQP2RK and blocked BLAS 3 code
682* in ZLAQP3RK.
683* 1) ZGEQP3RK, ZLAQP2RK, ZLAQP3RK: 2*N to store full and
684* partial column 2-norms.
685* 2) ZLAQP2RK: N+NRHS-1 to use in WORK array that is used
686* in ZLARF1F subroutine to apply an elementary reflector
687* from the left.
688* 3) ZLAQP3RK: NB*(N+NRHS) to use in the work array F that
689* is used to apply a block reflector from
690* the left.
691* 4) ZLAQP3RK: NB to use in the auxilixary array AUX.
692* Sizes (2) and ((3) + (4)) should intersect, therefore
693* TOTAL_WORK_SIZE = 2*N + NB*( N+NRHS+1 ), given NBMIN=2.
694*
695 lwkopt = 2*n + nb*( n+nrhs+1 )
696 END IF
697 work( 1 ) = dcmplx( lwkopt )
698*
699 IF( ( lwork.LT.iws ) .AND. .NOT.lquery ) THEN
700 info = -15
701 END IF
702 END IF
703*
704* NOTE: The optimal workspace size is returned in WORK(1), if
705* the input parameters M, N, NRHS, KMAX, LDA are valid.
706*
707 IF( info.NE.0 ) THEN
708 CALL xerbla( 'ZGEQP3RK', -info )
709 RETURN
710 ELSE IF( lquery ) THEN
711 RETURN
712 END IF
713*
714* Quick return if possible for M=0 or N=0.
715*
716 IF( minmn.EQ.0 ) THEN
717 k = 0
718 maxc2nrmk = zero
719 relmaxc2nrmk = zero
720 work( 1 ) = dcmplx( lwkopt )
721 RETURN
722 END IF
723*
724* ==================================================================
725*
726* Initialize column pivot array JPIV.
727*
728 DO j = 1, n
729 jpiv( j ) = j
730 END DO
731*
732* ==================================================================
733*
734* Initialize storage for partial and exact column 2-norms.
735* a) The elements WORK(1:N) are used to store partial column
736* 2-norms of the matrix A, and may decrease in each computation
737* step; initialize to the values of complete columns 2-norms.
738* b) The elements WORK(N+1:2*N) are used to store complete column
739* 2-norms of the matrix A, they are not changed during the
740* computation; initialize the values of complete columns 2-norms.
741*
742 DO j = 1, n
743 rwork( j ) = dznrm2( m, a( 1, j ), 1 )
744 rwork( n+j ) = rwork( j )
745 END DO
746*
747* ==================================================================
748*
749* Compute the pivot column index and the maximum column 2-norm
750* for the whole original matrix stored in A(1:M,1:N).
751*
752 kp1 = idamax( n, rwork( 1 ), 1 )
753*
754* ==================================================================.
755*
756 IF( disnan( maxc2nrm ) ) THEN
757*
758* Check if the matrix A contains NaN, set INFO parameter
759* to the column number where the first NaN is found and return
760* from the routine.
761*
762 k = 0
763 info = kp1
764*
765* Set MAXC2NRMK and RELMAXC2NRMK to NaN.
766*
767 maxc2nrmk = maxc2nrm
768 relmaxc2nrmk = maxc2nrm
769*
770* Array TAU is not set and contains undefined elements.
771*
772 work( 1 ) = dcmplx( lwkopt )
773 RETURN
774 END IF
775*
776* ===================================================================
777*
778 IF( maxc2nrm.EQ.zero ) THEN
779*
780* Check is the matrix A is a zero matrix, set array TAU and
781* return from the routine.
782*
783 k = 0
784 maxc2nrmk = zero
785 relmaxc2nrmk = zero
786*
787 DO j = 1, minmn
788 tau( j ) = czero
789 END DO
790*
791 work( 1 ) = dcmplx( lwkopt )
792 RETURN
793*
794 END IF
795*
796* ===================================================================
797*
798 hugeval = dlamch( 'Overflow' )
799*
800 IF( maxc2nrm.GT.hugeval ) THEN
801*
802* Check if the matrix A contains +Inf or -Inf, set INFO parameter
803* to the column number, where the first +/-Inf is found plus N,
804* and continue the computation.
805*
806 info = n + kp1
807*
808 END IF
809*
810* ==================================================================
811*
812* Quick return if possible for the case when the first
813* stopping criterion is satisfied, i.e. KMAX = 0.
814*
815 IF( kmax.EQ.0 ) THEN
816 k = 0
817 maxc2nrmk = maxc2nrm
818 relmaxc2nrmk = one
819 DO j = 1, minmn
820 tau( j ) = czero
821 END DO
822 work( 1 ) = dcmplx( lwkopt )
823 RETURN
824 END IF
825*
826* ==================================================================
827*
828 eps = dlamch('Epsilon')
829*
830* Adjust ABSTOL
831*
832 IF( abstol.GE.zero ) THEN
833 safmin = dlamch('Safe minimum')
834 abstol = max( abstol, two*safmin )
835 END IF
836*
837* Adjust RELTOL
838*
839 IF( reltol.GE.zero ) THEN
840 reltol = max( reltol, eps )
841 END IF
842*
843* ===================================================================
844*
845* JMAX is the maximum index of the column to be factorized,
846* which is also limited by the first stopping criterion KMAX.
847*
848 jmax = min( kmax, minmn )
849*
850* ===================================================================
851*
852* Quick return if possible for the case when the second or third
853* stopping criterion for the whole original matrix is satified,
854* i.e. MAXC2NRM <= ABSTOL or RELMAXC2NRM <= RELTOL
855* (which is ONE <= RELTOL).
856*
857 IF( maxc2nrm.LE.abstol .OR. one.LE.reltol ) THEN
858*
859 k = 0
860 maxc2nrmk = maxc2nrm
861 relmaxc2nrmk = one
862*
863 DO j = 1, minmn
864 tau( j ) = czero
865 END DO
866*
867 work( 1 ) = dcmplx( lwkopt )
868 RETURN
869 END IF
870*
871* ==================================================================
872* Factorize columns
873* ==================================================================
874*
875* Determine the block size.
876*
877 nbmin = 2
878 nx = 0
879*
880 IF( ( nb.GT.1 ) .AND. ( nb.LT.minmn ) ) THEN
881*
882* Determine when to cross over from blocked to unblocked code.
883* (for N less than NX, unblocked code should be used).
884*
885 nx = max( 0, ilaenv( ixover, 'ZGEQP3RK', ' ', m, n, -1,
886 $ -1 ) )
887*
888 IF( nx.LT.minmn ) THEN
889*
890* Determine if workspace is large enough for blocked code.
891*
892 IF( lwork.LT.lwkopt ) THEN
893*
894* Not enough workspace to use optimal block size that
895* is currently stored in NB.
896* Reduce NB and determine the minimum value of NB.
897*
898 nb = ( lwork-2*n ) / ( n+1 )
899 nbmin = max( 2, ilaenv( inbmin, 'ZGEQP3RK', ' ', m, n,
900 $ -1, -1 ) )
901*
902 END IF
903 END IF
904 END IF
905*
906* ==================================================================
907*
908* DONE is the boolean flag to rerpresent the case when the
909* factorization completed in the block factorization routine,
910* before the end of the block.
911*
912 done = .false.
913*
914* J is the column index.
915*
916 j = 1
917*
918* (1) Use blocked code initially.
919*
920* JMAXB is the maximum column index of the block, when the
921* blocked code is used, is also limited by the first stopping
922* criterion KMAX.
923*
924 jmaxb = min( kmax, minmn - nx )
925*
926 IF( nb.GE.nbmin .AND. nb.LT.jmax .AND. jmaxb.GT.0 ) THEN
927*
928* Loop over the column blocks of the matrix A(1:M,1:JMAXB). Here:
929* J is the column index of a column block;
930* JB is the column block size to pass to block factorization
931* routine in a loop step;
932* JBF is the number of columns that were actually factorized
933* that was returned by the block factorization routine
934* in a loop step, JBF <= JB;
935* N_SUB is the number of columns in the submatrix;
936* IOFFSET is the number of rows that should not be factorized.
937*
938 DO WHILE( j.LE.jmaxb )
939*
940 jb = min( nb, jmaxb-j+1 )
941 n_sub = n-j+1
942 ioffset = j-1
943*
944* Factorize JB columns among the columns A(J:N).
945*
946 CALL zlaqp3rk( m, n_sub, nrhs, ioffset, jb, abstol,
947 $ reltol, kp1, maxc2nrm, a( 1, j ), lda,
948 $ done, jbf, maxc2nrmk, relmaxc2nrmk,
949 $ jpiv( j ), tau( j ),
950 $ rwork( j ), rwork( n+j ),
951 $ work( 1 ), work( jb+1 ),
952 $ n+nrhs-j+1, iwork, iinfo )
953*
954* Set INFO on the first occurence of Inf.
955*
956 IF( iinfo.GT.n_sub .AND. info.EQ.0 ) THEN
957 info = 2*ioffset + iinfo
958 END IF
959*
960 IF( done ) THEN
961*
962* Either the submatrix is zero before the end of the
963* column block, or ABSTOL or RELTOL criterion is
964* satisfied before the end of the column block, we can
965* return from the routine. Perform the following before
966* returning:
967* a) Set the number of factorized columns K,
968* K = IOFFSET + JBF from the last call of blocked
969* routine.
970* NOTE: 1) MAXC2NRMK and RELMAXC2NRMK are returned
971* by the block factorization routine;
972* 2) The remaining TAUs are set to ZERO by the
973* block factorization routine.
974*
975 k = ioffset + jbf
976*
977* Set INFO on the first occurrence of NaN, NaN takes
978* prcedence over Inf.
979*
980 IF( iinfo.LE.n_sub .AND. iinfo.GT.0 ) THEN
981 info = ioffset + iinfo
982 END IF
983*
984* Return from the routine.
985*
986 work( 1 ) = dcmplx( lwkopt )
987*
988 RETURN
989*
990 END IF
991*
992 j = j + jbf
993*
994 END DO
995*
996 END IF
997*
998* Use unblocked code to factor the last or only block.
999* J = JMAX+1 means we factorized the maximum possible number of
1000* columns, that is in ELSE clause we need to compute
1001* the MAXC2NORM and RELMAXC2NORM to return after we processed
1002* the blocks.
1003*
1004 IF( j.LE.jmax ) THEN
1005*
1006* N_SUB is the number of columns in the submatrix;
1007* IOFFSET is the number of rows that should not be factorized.
1008*
1009 n_sub = n-j+1
1010 ioffset = j-1
1011*
1012 CALL zlaqp2rk( m, n_sub, nrhs, ioffset, jmax-j+1,
1013 $ abstol, reltol, kp1, maxc2nrm, a( 1, j ), lda,
1014 $ kf, maxc2nrmk, relmaxc2nrmk, jpiv( j ),
1015 $ tau( j ), rwork( j ), rwork( n+j ),
1016 $ work( 1 ), iinfo )
1017*
1018* ABSTOL or RELTOL criterion is satisfied when the number of
1019* the factorized columns KF is smaller then the number
1020* of columns JMAX-J+1 supplied to be factorized by the
1021* unblocked routine, we can return from
1022* the routine. Perform the following before returning:
1023* a) Set the number of factorized columns K,
1024* b) MAXC2NRMK and RELMAXC2NRMK are returned by the
1025* unblocked factorization routine above.
1026*
1027 k = j - 1 + kf
1028*
1029* Set INFO on the first exception occurence.
1030*
1031* Set INFO on the first exception occurence of Inf or NaN,
1032* (NaN takes precedence over Inf).
1033*
1034 IF( iinfo.GT.n_sub .AND. info.EQ.0 ) THEN
1035 info = 2*ioffset + iinfo
1036 ELSE IF( iinfo.LE.n_sub .AND. iinfo.GT.0 ) THEN
1037 info = ioffset + iinfo
1038 END IF
1039*
1040 ELSE
1041*
1042* Compute the return values for blocked code.
1043*
1044* Set the number of factorized columns if the unblocked routine
1045* was not called.
1046*
1047 k = jmax
1048*
1049* If there exits a residual matrix after the blocked code:
1050* 1) compute the values of MAXC2NRMK, RELMAXC2NRMK of the
1051* residual matrix, otherwise set them to ZERO;
1052* 2) Set TAU(K+1:MINMN) to ZERO.
1053*
1054 IF( k.LT.minmn ) THEN
1055 jmaxc2nrm = k + idamax( n-k, rwork( k+1 ), 1 )
1056 maxc2nrmk = rwork( jmaxc2nrm )
1057 IF( k.EQ.0 ) THEN
1058 relmaxc2nrmk = one
1059 ELSE
1060 relmaxc2nrmk = maxc2nrmk / maxc2nrm
1061 END IF
1062*
1063 DO j = k + 1, minmn
1064 tau( j ) = czero
1065 END DO
1066*
1067 ELSE
1068 maxc2nrmk = zero
1069 relmaxc2nrmk = zero
1070*
1071 END IF
1072*
1073* END IF( J.LE.JMAX ) THEN
1074*
1075 END IF
1076*
1077 work( 1 ) = dcmplx( lwkopt )
1078*
1079 RETURN
1080*
1081* End of ZGEQP3RK
1082*
1083 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgeqp3rk(m, n, nrhs, kmax, abstol, reltol, a, lda, k, maxc2nrmk, relmaxc2nrmk, jpiv, tau, work, lwork, rwork, iwork, info)
ZGEQP3RK computes a truncated Householder QR factorization with column pivoting of a complex m-by-n m...
Definition zgeqp3rk.f:582
subroutine zlaqp2rk(m, n, nrhs, ioffset, kmax, abstol, reltol, kp1, maxc2nrm, a, lda, k, maxc2nrmk, relmaxc2nrmk, jpiv, tau, vn1, vn2, work, info)
ZLAQP2RK computes truncated QR factorization with column pivoting of a complex matrix block using Lev...
Definition zlaqp2rk.f:335
subroutine zlaqp3rk(m, n, nrhs, ioffset, nb, abstol, reltol, kp1, maxc2nrm, a, lda, done, kb, maxc2nrmk, relmaxc2nrmk, jpiv, tau, vn1, vn2, auxv, f, ldf, iwork, info)
ZLAQP3RK computes a step of truncated QR factorization with column pivoting of a complex m-by-n matri...
Definition zlaqp3rk.f:386