File:  [local] / rpl / lapack / lapack / dlahr2.f
Revision 1.19: download - view: text, annotated - select for diffs - revision graph
Mon Aug 7 08:38:54 2023 UTC (9 months ago) by bertrand
Branches: MAIN
CVS tags: rpl-4_1_35, rpl-4_1_34, HEAD
Première mise à jour de lapack et blas.

    1: *> \brief \b DLAHR2 reduces the specified number of first columns of a general rectangular matrix A so that elements below the specified subdiagonal are zero, and returns auxiliary matrices which are needed to apply the transformation to the unreduced part of A.
    2: *
    3: *  =========== DOCUMENTATION ===========
    4: *
    5: * Online html documentation available at
    6: *            http://www.netlib.org/lapack/explore-html/
    7: *
    8: *> \htmlonly
    9: *> Download DLAHR2 + dependencies
   10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlahr2.f">
   11: *> [TGZ]</a>
   12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlahr2.f">
   13: *> [ZIP]</a>
   14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlahr2.f">
   15: *> [TXT]</a>
   16: *> \endhtmlonly
   17: *
   18: *  Definition:
   19: *  ===========
   20: *
   21: *       SUBROUTINE DLAHR2( N, K, NB, A, LDA, TAU, T, LDT, Y, LDY )
   22: *
   23: *       .. Scalar Arguments ..
   24: *       INTEGER            K, LDA, LDT, LDY, N, NB
   25: *       ..
   26: *       .. Array Arguments ..
   27: *       DOUBLE PRECISION  A( LDA, * ), T( LDT, NB ), TAU( NB ),
   28: *      $                   Y( LDY, NB )
   29: *       ..
   30: *
   31: *
   32: *> \par Purpose:
   33: *  =============
   34: *>
   35: *> \verbatim
   36: *>
   37: *> DLAHR2 reduces the first NB columns of A real general n-BY-(n-k+1)
   38: *> matrix A so that elements below the k-th subdiagonal are zero. The
   39: *> reduction is performed by an orthogonal similarity transformation
   40: *> Q**T * A * Q. The routine returns the matrices V and T which determine
   41: *> Q as a block reflector I - V*T*V**T, and also the matrix Y = A * V * T.
   42: *>
   43: *> This is an auxiliary routine called by DGEHRD.
   44: *> \endverbatim
   45: *
   46: *  Arguments:
   47: *  ==========
   48: *
   49: *> \param[in] N
   50: *> \verbatim
   51: *>          N is INTEGER
   52: *>          The order of the matrix A.
   53: *> \endverbatim
   54: *>
   55: *> \param[in] K
   56: *> \verbatim
   57: *>          K is INTEGER
   58: *>          The offset for the reduction. Elements below the k-th
   59: *>          subdiagonal in the first NB columns are reduced to zero.
   60: *>          K < N.
   61: *> \endverbatim
   62: *>
   63: *> \param[in] NB
   64: *> \verbatim
   65: *>          NB is INTEGER
   66: *>          The number of columns to be reduced.
   67: *> \endverbatim
   68: *>
   69: *> \param[in,out] A
   70: *> \verbatim
   71: *>          A is DOUBLE PRECISION array, dimension (LDA,N-K+1)
   72: *>          On entry, the n-by-(n-k+1) general matrix A.
   73: *>          On exit, the elements on and above the k-th subdiagonal in
   74: *>          the first NB columns are overwritten with the corresponding
   75: *>          elements of the reduced matrix; the elements below the k-th
   76: *>          subdiagonal, with the array TAU, represent the matrix Q as a
   77: *>          product of elementary reflectors. The other columns of A are
   78: *>          unchanged. See Further Details.
   79: *> \endverbatim
   80: *>
   81: *> \param[in] LDA
   82: *> \verbatim
   83: *>          LDA is INTEGER
   84: *>          The leading dimension of the array A.  LDA >= max(1,N).
   85: *> \endverbatim
   86: *>
   87: *> \param[out] TAU
   88: *> \verbatim
   89: *>          TAU is DOUBLE PRECISION array, dimension (NB)
   90: *>          The scalar factors of the elementary reflectors. See Further
   91: *>          Details.
   92: *> \endverbatim
   93: *>
   94: *> \param[out] T
   95: *> \verbatim
   96: *>          T is DOUBLE PRECISION array, dimension (LDT,NB)
   97: *>          The upper triangular matrix T.
   98: *> \endverbatim
   99: *>
  100: *> \param[in] LDT
  101: *> \verbatim
  102: *>          LDT is INTEGER
  103: *>          The leading dimension of the array T.  LDT >= NB.
  104: *> \endverbatim
  105: *>
  106: *> \param[out] Y
  107: *> \verbatim
  108: *>          Y is DOUBLE PRECISION array, dimension (LDY,NB)
  109: *>          The n-by-nb matrix Y.
  110: *> \endverbatim
  111: *>
  112: *> \param[in] LDY
  113: *> \verbatim
  114: *>          LDY is INTEGER
  115: *>          The leading dimension of the array Y. LDY >= N.
  116: *> \endverbatim
  117: *
  118: *  Authors:
  119: *  ========
  120: *
  121: *> \author Univ. of Tennessee
  122: *> \author Univ. of California Berkeley
  123: *> \author Univ. of Colorado Denver
  124: *> \author NAG Ltd.
  125: *
  126: *> \ingroup doubleOTHERauxiliary
  127: *
  128: *> \par Further Details:
  129: *  =====================
  130: *>
  131: *> \verbatim
  132: *>
  133: *>  The matrix Q is represented as a product of nb elementary reflectors
  134: *>
  135: *>     Q = H(1) H(2) . . . H(nb).
  136: *>
  137: *>  Each H(i) has the form
  138: *>
  139: *>     H(i) = I - tau * v * v**T
  140: *>
  141: *>  where tau is a real scalar, and v is a real vector with
  142: *>  v(1:i+k-1) = 0, v(i+k) = 1; v(i+k+1:n) is stored on exit in
  143: *>  A(i+k+1:n,i), and tau in TAU(i).
  144: *>
  145: *>  The elements of the vectors v together form the (n-k+1)-by-nb matrix
  146: *>  V which is needed, with T and Y, to apply the transformation to the
  147: *>  unreduced part of the matrix, using an update of the form:
  148: *>  A := (I - V*T*V**T) * (A - Y*V**T).
  149: *>
  150: *>  The contents of A on exit are illustrated by the following example
  151: *>  with n = 7, k = 3 and nb = 2:
  152: *>
  153: *>     ( a   a   a   a   a )
  154: *>     ( a   a   a   a   a )
  155: *>     ( a   a   a   a   a )
  156: *>     ( h   h   a   a   a )
  157: *>     ( v1  h   a   a   a )
  158: *>     ( v1  v2  a   a   a )
  159: *>     ( v1  v2  a   a   a )
  160: *>
  161: *>  where a denotes an element of the original matrix A, h denotes a
  162: *>  modified element of the upper Hessenberg matrix H, and vi denotes an
  163: *>  element of the vector defining H(i).
  164: *>
  165: *>  This subroutine is a slight modification of LAPACK-3.0's DLAHRD
  166: *>  incorporating improvements proposed by Quintana-Orti and Van de
  167: *>  Gejin. Note that the entries of A(1:K,2:NB) differ from those
  168: *>  returned by the original LAPACK-3.0's DLAHRD routine. (This
  169: *>  subroutine is not backward compatible with LAPACK-3.0's DLAHRD.)
  170: *> \endverbatim
  171: *
  172: *> \par References:
  173: *  ================
  174: *>
  175: *>  Gregorio Quintana-Orti and Robert van de Geijn, "Improving the
  176: *>  performance of reduction to Hessenberg form," ACM Transactions on
  177: *>  Mathematical Software, 32(2):180-194, June 2006.
  178: *>
  179: *  =====================================================================
  180:       SUBROUTINE DLAHR2( N, K, NB, A, LDA, TAU, T, LDT, Y, LDY )
  181: *
  182: *  -- LAPACK auxiliary routine --
  183: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
  184: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  185: *
  186: *     .. Scalar Arguments ..
  187:       INTEGER            K, LDA, LDT, LDY, N, NB
  188: *     ..
  189: *     .. Array Arguments ..
  190:       DOUBLE PRECISION  A( LDA, * ), T( LDT, NB ), TAU( NB ),
  191:      $                   Y( LDY, NB )
  192: *     ..
  193: *
  194: *  =====================================================================
  195: *
  196: *     .. Parameters ..
  197:       DOUBLE PRECISION  ZERO, ONE
  198:       PARAMETER          ( ZERO = 0.0D+0,
  199:      $                     ONE = 1.0D+0 )
  200: *     ..
  201: *     .. Local Scalars ..
  202:       INTEGER            I
  203:       DOUBLE PRECISION  EI
  204: *     ..
  205: *     .. External Subroutines ..
  206:       EXTERNAL           DAXPY, DCOPY, DGEMM, DGEMV, DLACPY,
  207:      $                   DLARFG, DSCAL, DTRMM, DTRMV
  208: *     ..
  209: *     .. Intrinsic Functions ..
  210:       INTRINSIC          MIN
  211: *     ..
  212: *     .. Executable Statements ..
  213: *
  214: *     Quick return if possible
  215: *
  216:       IF( N.LE.1 )
  217:      $   RETURN
  218: *
  219:       DO 10 I = 1, NB
  220:          IF( I.GT.1 ) THEN
  221: *
  222: *           Update A(K+1:N,I)
  223: *
  224: *           Update I-th column of A - Y * V**T
  225: *
  226:             CALL DGEMV( 'NO TRANSPOSE', N-K, I-1, -ONE, Y(K+1,1), LDY,
  227:      $                  A( K+I-1, 1 ), LDA, ONE, A( K+1, I ), 1 )
  228: *
  229: *           Apply I - V * T**T * V**T to this column (call it b) from the
  230: *           left, using the last column of T as workspace
  231: *
  232: *           Let  V = ( V1 )   and   b = ( b1 )   (first I-1 rows)
  233: *                    ( V2 )             ( b2 )
  234: *
  235: *           where V1 is unit lower triangular
  236: *
  237: *           w := V1**T * b1
  238: *
  239:             CALL DCOPY( I-1, A( K+1, I ), 1, T( 1, NB ), 1 )
  240:             CALL DTRMV( 'Lower', 'Transpose', 'UNIT',
  241:      $                  I-1, A( K+1, 1 ),
  242:      $                  LDA, T( 1, NB ), 1 )
  243: *
  244: *           w := w + V2**T * b2
  245: *
  246:             CALL DGEMV( 'Transpose', N-K-I+1, I-1,
  247:      $                  ONE, A( K+I, 1 ),
  248:      $                  LDA, A( K+I, I ), 1, ONE, T( 1, NB ), 1 )
  249: *
  250: *           w := T**T * w
  251: *
  252:             CALL DTRMV( 'Upper', 'Transpose', 'NON-UNIT',
  253:      $                  I-1, T, LDT,
  254:      $                  T( 1, NB ), 1 )
  255: *
  256: *           b2 := b2 - V2*w
  257: *
  258:             CALL DGEMV( 'NO TRANSPOSE', N-K-I+1, I-1, -ONE,
  259:      $                  A( K+I, 1 ),
  260:      $                  LDA, T( 1, NB ), 1, ONE, A( K+I, I ), 1 )
  261: *
  262: *           b1 := b1 - V1*w
  263: *
  264:             CALL DTRMV( 'Lower', 'NO TRANSPOSE',
  265:      $                  'UNIT', I-1,
  266:      $                  A( K+1, 1 ), LDA, T( 1, NB ), 1 )
  267:             CALL DAXPY( I-1, -ONE, T( 1, NB ), 1, A( K+1, I ), 1 )
  268: *
  269:             A( K+I-1, I-1 ) = EI
  270:          END IF
  271: *
  272: *        Generate the elementary reflector H(I) to annihilate
  273: *        A(K+I+1:N,I)
  274: *
  275:          CALL DLARFG( N-K-I+1, A( K+I, I ), A( MIN( K+I+1, N ), I ), 1,
  276:      $                TAU( I ) )
  277:          EI = A( K+I, I )
  278:          A( K+I, I ) = ONE
  279: *
  280: *        Compute  Y(K+1:N,I)
  281: *
  282:          CALL DGEMV( 'NO TRANSPOSE', N-K, N-K-I+1,
  283:      $               ONE, A( K+1, I+1 ),
  284:      $               LDA, A( K+I, I ), 1, ZERO, Y( K+1, I ), 1 )
  285:          CALL DGEMV( 'Transpose', N-K-I+1, I-1,
  286:      $               ONE, A( K+I, 1 ), LDA,
  287:      $               A( K+I, I ), 1, ZERO, T( 1, I ), 1 )
  288:          CALL DGEMV( 'NO TRANSPOSE', N-K, I-1, -ONE,
  289:      $               Y( K+1, 1 ), LDY,
  290:      $               T( 1, I ), 1, ONE, Y( K+1, I ), 1 )
  291:          CALL DSCAL( N-K, TAU( I ), Y( K+1, I ), 1 )
  292: *
  293: *        Compute T(1:I,I)
  294: *
  295:          CALL DSCAL( I-1, -TAU( I ), T( 1, I ), 1 )
  296:          CALL DTRMV( 'Upper', 'No Transpose', 'NON-UNIT',
  297:      $               I-1, T, LDT,
  298:      $               T( 1, I ), 1 )
  299:          T( I, I ) = TAU( I )
  300: *
  301:    10 CONTINUE
  302:       A( K+NB, NB ) = EI
  303: *
  304: *     Compute Y(1:K,1:NB)
  305: *
  306:       CALL DLACPY( 'ALL', K, NB, A( 1, 2 ), LDA, Y, LDY )
  307:       CALL DTRMM( 'RIGHT', 'Lower', 'NO TRANSPOSE',
  308:      $            'UNIT', K, NB,
  309:      $            ONE, A( K+1, 1 ), LDA, Y, LDY )
  310:       IF( N.GT.K+NB )
  311:      $   CALL DGEMM( 'NO TRANSPOSE', 'NO TRANSPOSE', K,
  312:      $               NB, N-K-NB, ONE,
  313:      $               A( 1, 2+NB ), LDA, A( K+1+NB, 1 ), LDA, ONE, Y,
  314:      $               LDY )
  315:       CALL DTRMM( 'RIGHT', 'Upper', 'NO TRANSPOSE',
  316:      $            'NON-UNIT', K, NB,
  317:      $            ONE, T, LDT, Y, LDY )
  318: *
  319:       RETURN
  320: *
  321: *     End of DLAHR2
  322: *
  323:       END

CVSweb interface <joel.bertrand@systella.fr>