LAPACK  3.4.1
LAPACK: Linear Algebra PACKage
slaqps.f
Go to the documentation of this file.
00001 *> \brief \b SLAQPS
00002 *
00003 *  =========== DOCUMENTATION ===========
00004 *
00005 * Online html documentation available at 
00006 *            http://www.netlib.org/lapack/explore-html/ 
00007 *
00008 *> \htmlonly
00009 *> Download SLAQPS + dependencies 
00010 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/slaqps.f"> 
00011 *> [TGZ]</a> 
00012 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/slaqps.f"> 
00013 *> [ZIP]</a> 
00014 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/slaqps.f"> 
00015 *> [TXT]</a>
00016 *> \endhtmlonly 
00017 *
00018 *  Definition:
00019 *  ===========
00020 *
00021 *       SUBROUTINE SLAQPS( M, N, OFFSET, NB, KB, A, LDA, JPVT, TAU, VN1,
00022 *                          VN2, AUXV, F, LDF )
00023 * 
00024 *       .. Scalar Arguments ..
00025 *       INTEGER            KB, LDA, LDF, M, N, NB, OFFSET
00026 *       ..
00027 *       .. Array Arguments ..
00028 *       INTEGER            JPVT( * )
00029 *       REAL               A( LDA, * ), AUXV( * ), F( LDF, * ), TAU( * ),
00030 *      $                   VN1( * ), VN2( * )
00031 *       ..
00032 *  
00033 *
00034 *> \par Purpose:
00035 *  =============
00036 *>
00037 *> \verbatim
00038 *>
00039 *> SLAQPS computes a step of QR factorization with column pivoting
00040 *> of a real M-by-N matrix A by using Blas-3.  It tries to factorize
00041 *> NB columns from A starting from the row OFFSET+1, and updates all
00042 *> of the matrix with Blas-3 xGEMM.
00043 *>
00044 *> In some cases, due to catastrophic cancellations, it cannot
00045 *> factorize NB columns.  Hence, the actual number of factorized
00046 *> columns is returned in KB.
00047 *>
00048 *> Block A(1:OFFSET,1:N) is accordingly pivoted, but not factorized.
00049 *> \endverbatim
00050 *
00051 *  Arguments:
00052 *  ==========
00053 *
00054 *> \param[in] M
00055 *> \verbatim
00056 *>          M is INTEGER
00057 *>          The number of rows of the matrix A. M >= 0.
00058 *> \endverbatim
00059 *>
00060 *> \param[in] N
00061 *> \verbatim
00062 *>          N is INTEGER
00063 *>          The number of columns of the matrix A. N >= 0
00064 *> \endverbatim
00065 *>
00066 *> \param[in] OFFSET
00067 *> \verbatim
00068 *>          OFFSET is INTEGER
00069 *>          The number of rows of A that have been factorized in
00070 *>          previous steps.
00071 *> \endverbatim
00072 *>
00073 *> \param[in] NB
00074 *> \verbatim
00075 *>          NB is INTEGER
00076 *>          The number of columns to factorize.
00077 *> \endverbatim
00078 *>
00079 *> \param[out] KB
00080 *> \verbatim
00081 *>          KB is INTEGER
00082 *>          The number of columns actually factorized.
00083 *> \endverbatim
00084 *>
00085 *> \param[in,out] A
00086 *> \verbatim
00087 *>          A is REAL array, dimension (LDA,N)
00088 *>          On entry, the M-by-N matrix A.
00089 *>          On exit, block A(OFFSET+1:M,1:KB) is the triangular
00090 *>          factor obtained and block A(1:OFFSET,1:N) has been
00091 *>          accordingly pivoted, but no factorized.
00092 *>          The rest of the matrix, block A(OFFSET+1:M,KB+1:N) has
00093 *>          been updated.
00094 *> \endverbatim
00095 *>
00096 *> \param[in] LDA
00097 *> \verbatim
00098 *>          LDA is INTEGER
00099 *>          The leading dimension of the array A. LDA >= max(1,M).
00100 *> \endverbatim
00101 *>
00102 *> \param[in,out] JPVT
00103 *> \verbatim
00104 *>          JPVT is INTEGER array, dimension (N)
00105 *>          JPVT(I) = K <==> Column K of the full matrix A has been
00106 *>          permuted into position I in AP.
00107 *> \endverbatim
00108 *>
00109 *> \param[out] TAU
00110 *> \verbatim
00111 *>          TAU is REAL array, dimension (KB)
00112 *>          The scalar factors of the elementary reflectors.
00113 *> \endverbatim
00114 *>
00115 *> \param[in,out] VN1
00116 *> \verbatim
00117 *>          VN1 is REAL array, dimension (N)
00118 *>          The vector with the partial column norms.
00119 *> \endverbatim
00120 *>
00121 *> \param[in,out] VN2
00122 *> \verbatim
00123 *>          VN2 is REAL array, dimension (N)
00124 *>          The vector with the exact column norms.
00125 *> \endverbatim
00126 *>
00127 *> \param[in,out] AUXV
00128 *> \verbatim
00129 *>          AUXV is REAL array, dimension (NB)
00130 *>          Auxiliar vector.
00131 *> \endverbatim
00132 *>
00133 *> \param[in,out] F
00134 *> \verbatim
00135 *>          F is REAL array, dimension (LDF,NB)
00136 *>          Matrix F**T = L*Y**T*A.
00137 *> \endverbatim
00138 *>
00139 *> \param[in] LDF
00140 *> \verbatim
00141 *>          LDF is INTEGER
00142 *>          The leading dimension of the array F. LDF >= max(1,N).
00143 *> \endverbatim
00144 *
00145 *  Authors:
00146 *  ========
00147 *
00148 *> \author Univ. of Tennessee 
00149 *> \author Univ. of California Berkeley 
00150 *> \author Univ. of Colorado Denver 
00151 *> \author NAG Ltd. 
00152 *
00153 *> \date November 2011
00154 *
00155 *> \ingroup realOTHERauxiliary
00156 *
00157 *> \par Contributors:
00158 *  ==================
00159 *>
00160 *>    G. Quintana-Orti, Depto. de Informatica, Universidad Jaime I, Spain
00161 *>    X. Sun, Computer Science Dept., Duke University, USA
00162 *>
00163 *> \n
00164 *>  Partial column norm updating strategy modified on April 2011
00165 *>    Z. Drmac and Z. Bujanovic, Dept. of Mathematics,
00166 *>    University of Zagreb, Croatia.
00167 *
00168 *> \par References:
00169 *  ================
00170 *>
00171 *> LAPACK Working Note 176
00172 *
00173 *> \htmlonly
00174 *> <a href="http://www.netlib.org/lapack/lawnspdf/lawn176.pdf">[PDF]</a> 
00175 *> \endhtmlonly 
00176 *
00177 *  =====================================================================
00178       SUBROUTINE SLAQPS( M, N, OFFSET, NB, KB, A, LDA, JPVT, TAU, VN1,
00179      $                   VN2, AUXV, F, LDF )
00180 *
00181 *  -- LAPACK auxiliary routine (version 3.4.0) --
00182 *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
00183 *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
00184 *     November 2011
00185 *
00186 *     .. Scalar Arguments ..
00187       INTEGER            KB, LDA, LDF, M, N, NB, OFFSET
00188 *     ..
00189 *     .. Array Arguments ..
00190       INTEGER            JPVT( * )
00191       REAL               A( LDA, * ), AUXV( * ), F( LDF, * ), TAU( * ),
00192      $                   VN1( * ), VN2( * )
00193 *     ..
00194 *
00195 *  =====================================================================
00196 *
00197 *     .. Parameters ..
00198       REAL               ZERO, ONE
00199       PARAMETER          ( ZERO = 0.0E+0, ONE = 1.0E+0 )
00200 *     ..
00201 *     .. Local Scalars ..
00202       INTEGER            ITEMP, J, K, LASTRK, LSTICC, PVT, RK
00203       REAL               AKK, TEMP, TEMP2, TOL3Z
00204 *     ..
00205 *     .. External Subroutines ..
00206       EXTERNAL           SGEMM, SGEMV, SLARFG, SSWAP
00207 *     ..
00208 *     .. Intrinsic Functions ..
00209       INTRINSIC          ABS, MAX, MIN, NINT, REAL, SQRT
00210 *     ..
00211 *     .. External Functions ..
00212       INTEGER            ISAMAX
00213       REAL               SLAMCH, SNRM2
00214       EXTERNAL           ISAMAX, SLAMCH, SNRM2
00215 *     ..
00216 *     .. Executable Statements ..
00217 *
00218       LASTRK = MIN( M, N+OFFSET )
00219       LSTICC = 0
00220       K = 0
00221       TOL3Z = SQRT(SLAMCH('Epsilon'))
00222 *
00223 *     Beginning of while loop.
00224 *
00225    10 CONTINUE
00226       IF( ( K.LT.NB ) .AND. ( LSTICC.EQ.0 ) ) THEN
00227          K = K + 1
00228          RK = OFFSET + K
00229 *
00230 *        Determine ith pivot column and swap if necessary
00231 *
00232          PVT = ( K-1 ) + ISAMAX( N-K+1, VN1( K ), 1 )
00233          IF( PVT.NE.K ) THEN
00234             CALL SSWAP( M, A( 1, PVT ), 1, A( 1, K ), 1 )
00235             CALL SSWAP( K-1, F( PVT, 1 ), LDF, F( K, 1 ), LDF )
00236             ITEMP = JPVT( PVT )
00237             JPVT( PVT ) = JPVT( K )
00238             JPVT( K ) = ITEMP
00239             VN1( PVT ) = VN1( K )
00240             VN2( PVT ) = VN2( K )
00241          END IF
00242 *
00243 *        Apply previous Householder reflectors to column K:
00244 *        A(RK:M,K) := A(RK:M,K) - A(RK:M,1:K-1)*F(K,1:K-1)**T.
00245 *
00246          IF( K.GT.1 ) THEN
00247             CALL SGEMV( 'No transpose', M-RK+1, K-1, -ONE, A( RK, 1 ),
00248      $                  LDA, F( K, 1 ), LDF, ONE, A( RK, K ), 1 )
00249          END IF
00250 *
00251 *        Generate elementary reflector H(k).
00252 *
00253          IF( RK.LT.M ) THEN
00254             CALL SLARFG( M-RK+1, A( RK, K ), A( RK+1, K ), 1, TAU( K ) )
00255          ELSE
00256             CALL SLARFG( 1, A( RK, K ), A( RK, K ), 1, TAU( K ) )
00257          END IF
00258 *
00259          AKK = A( RK, K )
00260          A( RK, K ) = ONE
00261 *
00262 *        Compute Kth column of F:
00263 *
00264 *        Compute  F(K+1:N,K) := tau(K)*A(RK:M,K+1:N)**T*A(RK:M,K).
00265 *
00266          IF( K.LT.N ) THEN
00267             CALL SGEMV( 'Transpose', M-RK+1, N-K, TAU( K ),
00268      $                  A( RK, K+1 ), LDA, A( RK, K ), 1, ZERO,
00269      $                  F( K+1, K ), 1 )
00270          END IF
00271 *
00272 *        Padding F(1:K,K) with zeros.
00273 *
00274          DO 20 J = 1, K
00275             F( J, K ) = ZERO
00276    20    CONTINUE
00277 *
00278 *        Incremental updating of F:
00279 *        F(1:N,K) := F(1:N,K) - tau(K)*F(1:N,1:K-1)*A(RK:M,1:K-1)**T
00280 *                    *A(RK:M,K).
00281 *
00282          IF( K.GT.1 ) THEN
00283             CALL SGEMV( 'Transpose', M-RK+1, K-1, -TAU( K ), A( RK, 1 ),
00284      $                  LDA, A( RK, K ), 1, ZERO, AUXV( 1 ), 1 )
00285 *
00286             CALL SGEMV( 'No transpose', N, K-1, ONE, F( 1, 1 ), LDF,
00287      $                  AUXV( 1 ), 1, ONE, F( 1, K ), 1 )
00288          END IF
00289 *
00290 *        Update the current row of A:
00291 *        A(RK,K+1:N) := A(RK,K+1:N) - A(RK,1:K)*F(K+1:N,1:K)**T.
00292 *
00293          IF( K.LT.N ) THEN
00294             CALL SGEMV( 'No transpose', N-K, K, -ONE, F( K+1, 1 ), LDF,
00295      $                  A( RK, 1 ), LDA, ONE, A( RK, K+1 ), LDA )
00296          END IF
00297 *
00298 *        Update partial column norms.
00299 *
00300          IF( RK.LT.LASTRK ) THEN
00301             DO 30 J = K + 1, N
00302                IF( VN1( J ).NE.ZERO ) THEN
00303 *
00304 *                 NOTE: The following 4 lines follow from the analysis in
00305 *                 Lapack Working Note 176.
00306 *
00307                   TEMP = ABS( A( RK, J ) ) / VN1( J )
00308                   TEMP = MAX( ZERO, ( ONE+TEMP )*( ONE-TEMP ) )
00309                   TEMP2 = TEMP*( VN1( J ) / VN2( J ) )**2
00310                   IF( TEMP2 .LE. TOL3Z ) THEN
00311                      VN2( J ) = REAL( LSTICC )
00312                      LSTICC = J
00313                   ELSE
00314                      VN1( J ) = VN1( J )*SQRT( TEMP )
00315                   END IF
00316                END IF
00317    30       CONTINUE
00318          END IF
00319 *
00320          A( RK, K ) = AKK
00321 *
00322 *        End of while loop.
00323 *
00324          GO TO 10
00325       END IF
00326       KB = K
00327       RK = OFFSET + KB
00328 *
00329 *     Apply the block reflector to the rest of the matrix:
00330 *     A(OFFSET+KB+1:M,KB+1:N) := A(OFFSET+KB+1:M,KB+1:N) -
00331 *                         A(OFFSET+KB+1:M,1:KB)*F(KB+1:N,1:KB)**T.
00332 *
00333       IF( KB.LT.MIN( N, M-OFFSET ) ) THEN
00334          CALL SGEMM( 'No transpose', 'Transpose', M-RK, N-KB, KB, -ONE,
00335      $               A( RK+1, 1 ), LDA, F( KB+1, 1 ), LDF, ONE,
00336      $               A( RK+1, KB+1 ), LDA )
00337       END IF
00338 *
00339 *     Recomputation of difficult columns.
00340 *
00341    40 CONTINUE
00342       IF( LSTICC.GT.0 ) THEN
00343          ITEMP = NINT( VN2( LSTICC ) )
00344          VN1( LSTICC ) = SNRM2( M-RK, A( RK+1, LSTICC ), 1 )
00345 *
00346 *        NOTE: The computation of VN1( LSTICC ) relies on the fact that 
00347 *        SNRM2 does not fail on vectors with norm below the value of
00348 *        SQRT(DLAMCH('S')) 
00349 *
00350          VN2( LSTICC ) = VN1( LSTICC )
00351          LSTICC = ITEMP
00352          GO TO 40
00353       END IF
00354 *
00355       RETURN
00356 *
00357 *     End of SLAQPS
00358 *
00359       END
 All Files Functions