LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
slaqp2rk.f
Go to the documentation of this file.
1*> \brief \b SLAQP2RK computes truncated QR factorization with column pivoting of a real matrix block using Level 2 BLAS and overwrites a real m-by-nrhs matrix B with Q**T * B.
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download SLAQP2RK + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/slaqp2rk.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/slaqp2rk.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/slaqp2rk.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE SLAQP2RK( M, N, NRHS, IOFFSET, KMAX, ABSTOL, RELTOL,
20* $ KP1, MAXC2NRM, A, LDA, K, MAXC2NRMK,
21* $ RELMAXC2NRMK, JPIV, TAU, VN1, VN2, WORK,
22* $ INFO )
23* IMPLICIT NONE
24*
25* .. Scalar Arguments ..
26* INTEGER INFO, IOFFSET, KP1, K, KMAX, LDA, M, N, NRHS
27* REAL ABSTOL, MAXC2NRM, MAXC2NRMK, RELMAXC2NRMK,
28* $ RELTOL
29* ..
30* .. Array Arguments ..
31* INTEGER JPIV( * )
32* REAL A( LDA, * ), TAU( * ), VN1( * ), VN2( * ),
33* $ WORK( * )
34* ..
35*
36*
37*> \par Purpose:
38* =============
39*>
40*> \verbatim
41*>
42*> SLAQP2RK computes a truncated (rank K) or full rank Householder QR
43*> factorization with column pivoting of a real matrix
44*> block A(IOFFSET+1:M,1:N) as
45*>
46*> A * P(K) = Q(K) * R(K).
47*>
48*> The routine uses Level 2 BLAS. The block A(1:IOFFSET,1:N)
49*> is accordingly pivoted, but not factorized.
50*>
51*> The routine also overwrites the right-hand-sides matrix block B
52*> stored in A(IOFFSET+1:M,N+1:N+NRHS) with Q(K)**T * B.
53*> \endverbatim
54*
55* Arguments:
56* ==========
57*
58*> \param[in] M
59*> \verbatim
60*> M is INTEGER
61*> The number of rows of the matrix A. M >= 0.
62*> \endverbatim
63*>
64*> \param[in] N
65*> \verbatim
66*> N is INTEGER
67*> The number of columns of the matrix A. N >= 0.
68*> \endverbatim
69*>
70*> \param[in] NRHS
71*> \verbatim
72*> NRHS is INTEGER
73*> The number of right hand sides, i.e., the number of
74*> columns of the matrix B. NRHS >= 0.
75*> \endverbatim
76*>
77*> \param[in] IOFFSET
78*> \verbatim
79*> IOFFSET is INTEGER
80*> The number of rows of the matrix A that must be pivoted
81*> but not factorized. IOFFSET >= 0.
82*>
83*> IOFFSET also represents the number of columns of the whole
84*> original matrix A_orig that have been factorized
85*> in the previous steps.
86*> \endverbatim
87*>
88*> \param[in] KMAX
89*> \verbatim
90*> KMAX is INTEGER
91*>
92*> The first factorization stopping criterion. KMAX >= 0.
93*>
94*> The maximum number of columns of the matrix A to factorize,
95*> i.e. the maximum factorization rank.
96*>
97*> a) If KMAX >= min(M-IOFFSET,N), then this stopping
98*> criterion is not used, factorize columns
99*> depending on ABSTOL and RELTOL.
100*>
101*> b) If KMAX = 0, then this stopping criterion is
102*> satisfied on input and the routine exits immediately.
103*> This means that the factorization is not performed,
104*> the matrices A and B and the arrays TAU, IPIV
105*> are not modified.
106*> \endverbatim
107*>
108*> \param[in] ABSTOL
109*> \verbatim
110*> ABSTOL is DOUBLE PRECISION, cannot be NaN.
111*>
112*> The second factorization stopping criterion.
113*>
114*> The absolute tolerance (stopping threshold) for
115*> maximum column 2-norm of the residual matrix.
116*> The algorithm converges (stops the factorization) when
117*> the maximum column 2-norm of the residual matrix
118*> is less than or equal to ABSTOL.
119*>
120*> a) If ABSTOL < 0.0, then this stopping criterion is not
121*> used, the routine factorizes columns depending
122*> on KMAX and RELTOL.
123*> This includes the case ABSTOL = -Inf.
124*>
125*> b) If 0.0 <= ABSTOL then the input value
126*> of ABSTOL is used.
127*> \endverbatim
128*>
129*> \param[in] RELTOL
130*> \verbatim
131*> RELTOL is DOUBLE PRECISION, cannot be NaN.
132*>
133*> The third factorization stopping criterion.
134*>
135*> The tolerance (stopping threshold) for the ratio of the
136*> maximum column 2-norm of the residual matrix to the maximum
137*> column 2-norm of the original matrix A_orig. The algorithm
138*> converges (stops the factorization), when this ratio is
139*> less than or equal to RELTOL.
140*>
141*> a) If RELTOL < 0.0, then this stopping criterion is not
142*> used, the routine factorizes columns depending
143*> on KMAX and ABSTOL.
144*> This includes the case RELTOL = -Inf.
145*>
146*> d) If 0.0 <= RELTOL then the input value of RELTOL
147*> is used.
148*> \endverbatim
149*>
150*> \param[in] KP1
151*> \verbatim
152*> KP1 is INTEGER
153*> The index of the column with the maximum 2-norm in
154*> the whole original matrix A_orig determined in the
155*> main routine SGEQP3RK. 1 <= KP1 <= N_orig_mat.
156*> \endverbatim
157*>
158*> \param[in] MAXC2NRM
159*> \verbatim
160*> MAXC2NRM is DOUBLE PRECISION
161*> The maximum column 2-norm of the whole original
162*> matrix A_orig computed in the main routine SGEQP3RK.
163*> MAXC2NRM >= 0.
164*> \endverbatim
165*>
166*> \param[in,out] A
167*> \verbatim
168*> A is REAL array, dimension (LDA,N+NRHS)
169*> On entry:
170*> the M-by-N matrix A and M-by-NRHS matrix B, as in
171*>
172*> N NRHS
173*> array_A = M [ mat_A, mat_B ]
174*>
175*> On exit:
176*> 1. The elements in block A(IOFFSET+1:M,1:K) below
177*> the diagonal together with the array TAU represent
178*> the orthogonal matrix Q(K) as a product of elementary
179*> reflectors.
180*> 2. The upper triangular block of the matrix A stored
181*> in A(IOFFSET+1:M,1:K) is the triangular factor obtained.
182*> 3. The block of the matrix A stored in A(1:IOFFSET,1:N)
183*> has been accordingly pivoted, but not factorized.
184*> 4. The rest of the array A, block A(IOFFSET+1:M,K+1:N+NRHS).
185*> The left part A(IOFFSET+1:M,K+1:N) of this block
186*> contains the residual of the matrix A, and,
187*> if NRHS > 0, the right part of the block
188*> A(IOFFSET+1:M,N+1:N+NRHS) contains the block of
189*> the right-hand-side matrix B. Both these blocks have been
190*> updated by multiplication from the left by Q(K)**T.
191*> \endverbatim
192*>
193*> \param[in] LDA
194*> \verbatim
195*> LDA is INTEGER
196*> The leading dimension of the array A. LDA >= max(1,M).
197*> \endverbatim
198*>
199*> \param[out] K
200*> \verbatim
201*> K is INTEGER
202*> Factorization rank of the matrix A, i.e. the rank of
203*> the factor R, which is the same as the number of non-zero
204*> rows of the factor R. 0 <= K <= min(M-IOFFSET,KMAX,N).
205*>
206*> K also represents the number of non-zero Householder
207*> vectors.
208*> \endverbatim
209*>
210*> \param[out] MAXC2NRMK
211*> \verbatim
212*> MAXC2NRMK is DOUBLE PRECISION
213*> The maximum column 2-norm of the residual matrix,
214*> when the factorization stopped at rank K. MAXC2NRMK >= 0.
215*> \endverbatim
216*>
217*> \param[out] RELMAXC2NRMK
218*> \verbatim
219*> RELMAXC2NRMK is DOUBLE PRECISION
220*> The ratio MAXC2NRMK / MAXC2NRM of the maximum column
221*> 2-norm of the residual matrix (when the factorization
222*> stopped at rank K) to the maximum column 2-norm of the
223*> whole original matrix A. RELMAXC2NRMK >= 0.
224*> \endverbatim
225*>
226*> \param[out] JPIV
227*> \verbatim
228*> JPIV is INTEGER array, dimension (N)
229*> Column pivot indices, for 1 <= j <= N, column j
230*> of the matrix A was interchanged with column JPIV(j).
231*> \endverbatim
232*>
233*> \param[out] TAU
234*> \verbatim
235*> TAU is REAL array, dimension (min(M-IOFFSET,N))
236*> The scalar factors of the elementary reflectors.
237*> \endverbatim
238*>
239*> \param[in,out] VN1
240*> \verbatim
241*> VN1 is REAL array, dimension (N)
242*> The vector with the partial column norms.
243*> \endverbatim
244*>
245*> \param[in,out] VN2
246*> \verbatim
247*> VN2 is REAL array, dimension (N)
248*> The vector with the exact column norms.
249*> \endverbatim
250*>
251*> \param[out] WORK
252*> \verbatim
253*> WORK is REAL array, dimension (N-1)
254*> Used in SLARF1F subroutine to apply an elementary
255*> reflector from the left.
256*> \endverbatim
257*>
258*> \param[out] INFO
259*> \verbatim
260*> INFO is INTEGER
261*> 1) INFO = 0: successful exit.
262*> 2) If INFO = j_1, where 1 <= j_1 <= N, then NaN was
263*> detected and the routine stops the computation.
264*> The j_1-th column of the matrix A or the j_1-th
265*> element of array TAU contains the first occurrence
266*> of NaN in the factorization step K+1 ( when K columns
267*> have been factorized ).
268*>
269*> On exit:
270*> K is set to the number of
271*> factorized columns without
272*> exception.
273*> MAXC2NRMK is set to NaN.
274*> RELMAXC2NRMK is set to NaN.
275*> TAU(K+1:min(M,N)) is not set and contains undefined
276*> elements. If j_1=K+1, TAU(K+1)
277*> may contain NaN.
278*> 3) If INFO = j_2, where N+1 <= j_2 <= 2*N, then no NaN
279*> was detected, but +Inf (or -Inf) was detected and
280*> the routine continues the computation until completion.
281*> The (j_2-N)-th column of the matrix A contains the first
282*> occurrence of +Inf (or -Inf) in the factorization
283*> step K+1 ( when K columns have been factorized ).
284*> \endverbatim
285*
286* Authors:
287* ========
288*
289*> \author Univ. of Tennessee
290*> \author Univ. of California Berkeley
291*> \author Univ. of Colorado Denver
292*> \author NAG Ltd.
293*
294*> \ingroup laqp2rk
295*
296*> \par References:
297* ================
298*> [1] A Level 3 BLAS QR factorization algorithm with column pivoting developed in 1996.
299*> G. Quintana-Orti, Depto. de Informatica, Universidad Jaime I, Spain.
300*> X. Sun, Computer Science Dept., Duke University, USA.
301*> C. H. Bischof, Math. and Comp. Sci. Div., Argonne National Lab, USA.
302*> A BLAS-3 version of the QR factorization with column pivoting.
303*> LAPACK Working Note 114
304*> <a href="https://www.netlib.org/lapack/lawnspdf/lawn114.pdf">https://www.netlib.org/lapack/lawnspdf/lawn114.pdf</a>
305*> and in
306*> SIAM J. Sci. Comput., 19(5):1486-1494, Sept. 1998.
307*> <a href="https://doi.org/10.1137/S1064827595296732">https://doi.org/10.1137/S1064827595296732</a>
308*>
309*> [2] A partial column norm updating strategy developed in 2006.
310*> Z. Drmac and Z. Bujanovic, Dept. of Math., University of Zagreb, Croatia.
311*> On the failure of rank revealing QR factorization software – a case study.
312*> LAPACK Working Note 176.
313*> <a href="http://www.netlib.org/lapack/lawnspdf/lawn176.pdf">http://www.netlib.org/lapack/lawnspdf/lawn176.pdf</a>
314*> and in
315*> ACM Trans. Math. Softw. 35, 2, Article 12 (July 2008), 28 pages.
316*> <a href="https://doi.org/10.1145/1377612.1377616">https://doi.org/10.1145/1377612.1377616</a>
317*
318*> \par Contributors:
319* ==================
320*>
321*> \verbatim
322*>
323*> November 2023, Igor Kozachenko, James Demmel,
324*> EECS Department,
325*> University of California, Berkeley, USA.
326*>
327*> \endverbatim
328*
329* =====================================================================
330 SUBROUTINE slaqp2rk( M, N, NRHS, IOFFSET, KMAX, ABSTOL, RELTOL,
331 $ KP1, MAXC2NRM, A, LDA, K, MAXC2NRMK,
332 $ RELMAXC2NRMK, JPIV, TAU, VN1, VN2, WORK,
333 $ INFO )
334 IMPLICIT NONE
335*
336* -- LAPACK auxiliary routine --
337* -- LAPACK is a software package provided by Univ. of Tennessee, --
338* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
339*
340* .. Scalar Arguments ..
341 INTEGER INFO, IOFFSET, KP1, K, KMAX, LDA, M, N, NRHS
342 REAL ABSTOL, MAXC2NRM, MAXC2NRMK, RELMAXC2NRMK,
343 $ RELTOL
344* ..
345* .. Array Arguments ..
346 INTEGER JPIV( * )
347 REAL A( LDA, * ), TAU( * ), VN1( * ), VN2( * ),
348 $ WORK( * )
349* ..
350*
351* =====================================================================
352*
353* .. Parameters ..
354 REAL ZERO, ONE
355 PARAMETER ( ZERO = 0.0e+0, one = 1.0e+0 )
356* ..
357* .. Local Scalars ..
358 INTEGER I, ITEMP, J, JMAXC2NRM, KK, KP, MINMNFACT,
359 $ MINMNUPDT
360 REAL HUGEVAL, TEMP, TEMP2, TOL3Z
361* ..
362* .. External Subroutines ..
363 EXTERNAL slarf1f, slarfg, sswap
364* ..
365* .. Intrinsic Functions ..
366 INTRINSIC abs, max, min, sqrt
367* ..
368* .. External Functions ..
369 LOGICAL SISNAN
370 INTEGER ISAMAX
371 REAL SLAMCH, SNRM2
372 EXTERNAL sisnan, slamch, isamax, snrm2
373* ..
374* .. Executable Statements ..
375*
376* Initialize INFO
377*
378 info = 0
379*
380* MINMNFACT in the smallest dimension of the submatrix
381* A(IOFFSET+1:M,1:N) to be factorized.
382*
383* MINMNUPDT is the smallest dimension
384* of the subarray A(IOFFSET+1:M,1:N+NRHS) to be udated, which
385* contains the submatrices A(IOFFSET+1:M,1:N) and
386* B(IOFFSET+1:M,1:NRHS) as column blocks.
387*
388 minmnfact = min( m-ioffset, n )
389 minmnupdt = min( m-ioffset, n+nrhs )
390 kmax = min( kmax, minmnfact )
391 tol3z = sqrt( slamch( 'Epsilon' ) )
392 hugeval = slamch( 'Overflow' )
393*
394* Compute the factorization, KK is the lomn loop index.
395*
396 DO kk = 1, kmax
397*
398 i = ioffset + kk
399*
400 IF( i.EQ.1 ) THEN
401*
402* ============================================================
403*
404* We are at the first column of the original whole matrix A,
405* therefore we use the computed KP1 and MAXC2NRM from the
406* main routine.
407*
408
409 kp = kp1
410*
411* ============================================================
412*
413 ELSE
414*
415* ============================================================
416*
417* Determine the pivot column in KK-th step, i.e. the index
418* of the column with the maximum 2-norm in the
419* submatrix A(I:M,K:N).
420*
421 kp = ( kk-1 ) + isamax( n-kk+1, vn1( kk ), 1 )
422*
423* Determine the maximum column 2-norm and the relative maximum
424* column 2-norm of the submatrix A(I:M,KK:N) in step KK.
425* RELMAXC2NRMK will be computed later, after somecondition
426* checks on MAXC2NRMK.
427*
428 maxc2nrmk = vn1( kp )
429*
430* ============================================================
431*
432* Check if the submatrix A(I:M,KK:N) contains NaN, and set
433* INFO parameter to the column number, where the first NaN
434* is found and return from the routine.
435* We need to check the condition only if the
436* column index (same as row index) of the original whole
437* matrix is larger than 1, since the condition for whole
438* original matrix is checked in the main routine.
439*
440 IF( sisnan( maxc2nrmk ) ) THEN
441*
442* Set K, the number of factorized columns.
443* that are not zero.
444*
445 k = kk - 1
446 info = k + kp
447*
448* Set RELMAXC2NRMK to NaN.
449*
450 relmaxc2nrmk = maxc2nrmk
451*
452* Array TAU(K+1:MINMNFACT) is not set and contains
453* undefined elements.
454*
455 RETURN
456 END IF
457*
458* ============================================================
459*
460* Quick return, if the submatrix A(I:M,KK:N) is
461* a zero matrix.
462* We need to check the condition only if the
463* column index (same as row index) of the original whole
464* matrix is larger than 1, since the condition for whole
465* original matrix is checked in the main routine.
466*
467 IF( maxc2nrmk.EQ.zero ) THEN
468*
469* Set K, the number of factorized columns.
470* that are not zero.
471*
472 k = kk - 1
473 relmaxc2nrmk = zero
474*
475* Set TAUs corresponding to the columns that were not
476* factorized to ZERO, i.e. set TAU(KK:MINMNFACT) to ZERO.
477*
478 DO j = kk, minmnfact
479 tau( j ) = zero
480 END DO
481*
482* Return from the routine.
483*
484 RETURN
485*
486 END IF
487*
488* ============================================================
489*
490* Check if the submatrix A(I:M,KK:N) contains Inf,
491* set INFO parameter to the column number, where
492* the first Inf is found plus N, and continue
493* the computation.
494* We need to check the condition only if the
495* column index (same as row index) of the original whole
496* matrix is larger than 1, since the condition for whole
497* original matrix is checked in the main routine.
498*
499 IF( info.EQ.0 .AND. maxc2nrmk.GT.hugeval ) THEN
500 info = n + kk - 1 + kp
501 END IF
502*
503* ============================================================
504*
505* Test for the second and third stopping criteria.
506* NOTE: There is no need to test for ABSTOL >= ZERO, since
507* MAXC2NRMK is non-negative. Similarly, there is no need
508* to test for RELTOL >= ZERO, since RELMAXC2NRMK is
509* non-negative.
510* We need to check the condition only if the
511* column index (same as row index) of the original whole
512* matrix is larger than 1, since the condition for whole
513* original matrix is checked in the main routine.
514
515 relmaxc2nrmk = maxc2nrmk / maxc2nrm
516*
517 IF( maxc2nrmk.LE.abstol .OR. relmaxc2nrmk.LE.reltol ) THEN
518*
519* Set K, the number of factorized columns.
520*
521 k = kk - 1
522*
523* Set TAUs corresponding to the columns that were not
524* factorized to ZERO, i.e. set TAU(KK:MINMNFACT) to ZERO.
525*
526 DO j = kk, minmnfact
527 tau( j ) = zero
528 END DO
529*
530* Return from the routine.
531*
532 RETURN
533*
534 END IF
535*
536* ============================================================
537*
538* End ELSE of IF(I.EQ.1)
539*
540 END IF
541*
542* ===============================================================
543*
544* If the pivot column is not the first column of the
545* subblock A(1:M,KK:N):
546* 1) swap the KK-th column and the KP-th pivot column
547* in A(1:M,1:N);
548* 2) copy the KK-th element into the KP-th element of the partial
549* and exact 2-norm vectors VN1 and VN2. ( Swap is not needed
550* for VN1 and VN2 since we use the element with the index
551* larger than KK in the next loop step.)
552* 3) Save the pivot interchange with the indices relative to the
553* the original matrix A, not the block A(1:M,1:N).
554*
555 IF( kp.NE.kk ) THEN
556 CALL sswap( m, a( 1, kp ), 1, a( 1, kk ), 1 )
557 vn1( kp ) = vn1( kk )
558 vn2( kp ) = vn2( kk )
559 itemp = jpiv( kp )
560 jpiv( kp ) = jpiv( kk )
561 jpiv( kk ) = itemp
562 END IF
563*
564* Generate elementary reflector H(KK) using the column A(I:M,KK),
565* if the column has more than one element, otherwise
566* the elementary reflector would be an identity matrix,
567* and TAU(KK) = ZERO.
568*
569 IF( i.LT.m ) THEN
570 CALL slarfg( m-i+1, a( i, kk ), a( i+1, kk ), 1,
571 $ tau( kk ) )
572 ELSE
573 tau( kk ) = zero
574 END IF
575*
576* Check if TAU(KK) contains NaN, set INFO parameter
577* to the column number where NaN is found and return from
578* the routine.
579* NOTE: There is no need to check TAU(KK) for Inf,
580* since SLARFG cannot produce TAU(KK) or Householder vector
581* below the diagonal containing Inf. Only BETA on the diagonal,
582* returned by SLARFG can contain Inf, which requires
583* TAU(KK) to contain NaN. Therefore, this case of generating Inf
584* by SLARFG is covered by checking TAU(KK) for NaN.
585*
586 IF( sisnan( tau(kk) ) ) THEN
587 k = kk - 1
588 info = kk
589*
590* Set MAXC2NRMK and RELMAXC2NRMK to NaN.
591*
592 maxc2nrmk = tau( kk )
593 relmaxc2nrmk = tau( kk )
594*
595* Array TAU(KK:MINMNFACT) is not set and contains
596* undefined elements, except the first element TAU(KK) = NaN.
597*
598 RETURN
599 END IF
600*
601* Apply H(KK)**T to A(I:M,KK+1:N+NRHS) from the left.
602* ( If M >= N, then at KK = N there is no residual matrix,
603* i.e. no columns of A to update, only columns of B.
604* If M < N, then at KK = M-IOFFSET, I = M and we have a
605* one-row residual matrix in A and the elementary
606* reflector is a unit matrix, TAU(KK) = ZERO, i.e. no update
607* is needed for the residual matrix in A and the
608* right-hand-side-matrix in B.
609* Therefore, we update only if
610* KK < MINMNUPDT = min(M-IOFFSET, N+NRHS)
611* condition is satisfied, not only KK < N+NRHS )
612*
613 IF( kk.LT.minmnupdt ) THEN
614 CALL slarf1f( 'Left', m-i+1, n+nrhs-kk, a( i, kk ), 1,
615 $ tau( kk ), a( i, kk+1 ), lda, work( 1 ) )
616 END IF
617*
618 IF( kk.LT.minmnfact ) THEN
619*
620* Update the partial column 2-norms for the residual matrix,
621* only if the residual matrix A(I+1:M,KK+1:N) exists, i.e.
622* when KK < min(M-IOFFSET, N).
623*
624 DO j = kk + 1, n
625 IF( vn1( j ).NE.zero ) THEN
626*
627* NOTE: The following lines follow from the analysis in
628* Lapack Working Note 176.
629*
630 temp = one - ( abs( a( i, j ) ) / vn1( j ) )**2
631 temp = max( temp, zero )
632 temp2 = temp*( vn1( j ) / vn2( j ) )**2
633 IF( temp2 .LE. tol3z ) THEN
634*
635* Compute the column 2-norm for the partial
636* column A(I+1:M,J) by explicitly computing it,
637* and store it in both partial 2-norm vector VN1
638* and exact column 2-norm vector VN2.
639*
640 vn1( j ) = snrm2( m-i, a( i+1, j ), 1 )
641 vn2( j ) = vn1( j )
642*
643 ELSE
644*
645* Update the column 2-norm for the partial
646* column A(I+1:M,J) by removing one
647* element A(I,J) and store it in partial
648* 2-norm vector VN1.
649*
650 vn1( j ) = vn1( j )*sqrt( temp )
651*
652 END IF
653 END IF
654 END DO
655*
656 END IF
657*
658* End factorization loop
659*
660 END DO
661*
662* If we reached this point, all colunms have been factorized,
663* i.e. no condition was triggered to exit the routine.
664* Set the number of factorized columns.
665*
666 k = kmax
667*
668* We reached the end of the loop, i.e. all KMAX columns were
669* factorized, we need to set MAXC2NRMK and RELMAXC2NRMK before
670* we return.
671*
672 IF( k.LT.minmnfact ) THEN
673*
674 jmaxc2nrm = k + isamax( n-k, vn1( k+1 ), 1 )
675 maxc2nrmk = vn1( jmaxc2nrm )
676*
677 IF( k.EQ.0 ) THEN
678 relmaxc2nrmk = one
679 ELSE
680 relmaxc2nrmk = maxc2nrmk / maxc2nrm
681 END IF
682*
683 ELSE
684 maxc2nrmk = zero
685 relmaxc2nrmk = zero
686 END IF
687*
688* We reached the end of the loop, i.e. all KMAX columns were
689* factorized, set TAUs corresponding to the columns that were
690* not factorized to ZERO, i.e. TAU(K+1:MINMNFACT) set to ZERO.
691*
692 DO j = k + 1, minmnfact
693 tau( j ) = zero
694 END DO
695*
696 RETURN
697*
698* End of SLAQP2RK
699*
700 END
subroutine slarfg(n, alpha, x, incx, tau)
SLARFG generates an elementary reflector (Householder matrix).
Definition slarfg.f:104
subroutine sswap(n, sx, incx, sy, incy)
SSWAP
Definition sswap.f:82
subroutine slaqp2rk(m, n, nrhs, ioffset, kmax, abstol, reltol, kp1, maxc2nrm, a, lda, k, maxc2nrmk, relmaxc2nrmk, jpiv, tau, vn1, vn2, work, info)
SLAQP2RK computes truncated QR factorization with column pivoting of a real matrix block using Level ...
Definition slaqp2rk.f:334
subroutine slarf1f(side, m, n, v, incv, tau, c, ldc, work)
SLARF1F applies an elementary reflector to a general rectangular
Definition slarf1f.f:123