Annotation of rpl/lapack/lapack/zlarfg.f, revision 1.20

1.12      bertrand    1: *> \brief \b ZLARFG generates an elementary reflector (Householder matrix).
1.9       bertrand    2: *
                      3: *  =========== DOCUMENTATION ===========
                      4: *
1.16      bertrand    5: * Online html documentation available at
                      6: *            http://www.netlib.org/lapack/explore-html/
1.9       bertrand    7: *
                      8: *> \htmlonly
1.16      bertrand    9: *> Download ZLARFG + dependencies
                     10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zlarfg.f">
                     11: *> [TGZ]</a>
                     12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zlarfg.f">
                     13: *> [ZIP]</a>
                     14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zlarfg.f">
1.9       bertrand   15: *> [TXT]</a>
1.16      bertrand   16: *> \endhtmlonly
1.9       bertrand   17: *
                     18: *  Definition:
                     19: *  ===========
                     20: *
                     21: *       SUBROUTINE ZLARFG( N, ALPHA, X, INCX, TAU )
1.16      bertrand   22: *
1.9       bertrand   23: *       .. Scalar Arguments ..
                     24: *       INTEGER            INCX, N
                     25: *       COMPLEX*16         ALPHA, TAU
                     26: *       ..
                     27: *       .. Array Arguments ..
                     28: *       COMPLEX*16         X( * )
                     29: *       ..
1.16      bertrand   30: *
1.9       bertrand   31: *
                     32: *> \par Purpose:
                     33: *  =============
                     34: *>
                     35: *> \verbatim
                     36: *>
                     37: *> ZLARFG generates a complex elementary reflector H of order n, such
                     38: *> that
                     39: *>
                     40: *>       H**H * ( alpha ) = ( beta ),   H**H * H = I.
                     41: *>              (   x   )   (   0  )
                     42: *>
                     43: *> where alpha and beta are scalars, with beta real, and x is an
                     44: *> (n-1)-element complex vector. H is represented in the form
                     45: *>
                     46: *>       H = I - tau * ( 1 ) * ( 1 v**H ) ,
                     47: *>                     ( v )
                     48: *>
                     49: *> where tau is a complex scalar and v is a complex (n-1)-element
                     50: *> vector. Note that H is not hermitian.
                     51: *>
                     52: *> If the elements of x are all zero and alpha is real, then tau = 0
                     53: *> and H is taken to be the unit matrix.
                     54: *>
                     55: *> Otherwise  1 <= real(tau) <= 2  and  abs(tau-1) <= 1 .
                     56: *> \endverbatim
                     57: *
                     58: *  Arguments:
                     59: *  ==========
                     60: *
                     61: *> \param[in] N
                     62: *> \verbatim
                     63: *>          N is INTEGER
                     64: *>          The order of the elementary reflector.
                     65: *> \endverbatim
                     66: *>
                     67: *> \param[in,out] ALPHA
                     68: *> \verbatim
                     69: *>          ALPHA is COMPLEX*16
                     70: *>          On entry, the value alpha.
                     71: *>          On exit, it is overwritten with the value beta.
                     72: *> \endverbatim
                     73: *>
                     74: *> \param[in,out] X
                     75: *> \verbatim
                     76: *>          X is COMPLEX*16 array, dimension
                     77: *>                         (1+(N-2)*abs(INCX))
                     78: *>          On entry, the vector x.
                     79: *>          On exit, it is overwritten with the vector v.
                     80: *> \endverbatim
                     81: *>
                     82: *> \param[in] INCX
                     83: *> \verbatim
                     84: *>          INCX is INTEGER
                     85: *>          The increment between elements of X. INCX > 0.
                     86: *> \endverbatim
                     87: *>
                     88: *> \param[out] TAU
                     89: *> \verbatim
                     90: *>          TAU is COMPLEX*16
                     91: *>          The value tau.
                     92: *> \endverbatim
                     93: *
                     94: *  Authors:
                     95: *  ========
                     96: *
1.16      bertrand   97: *> \author Univ. of Tennessee
                     98: *> \author Univ. of California Berkeley
                     99: *> \author Univ. of Colorado Denver
                    100: *> \author NAG Ltd.
1.9       bertrand  101: *
                    102: *> \ingroup complex16OTHERauxiliary
                    103: *
                    104: *  =====================================================================
1.1       bertrand  105:       SUBROUTINE ZLARFG( N, ALPHA, X, INCX, TAU )
                    106: *
1.20    ! bertrand  107: *  -- LAPACK auxiliary routine --
1.1       bertrand  108: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
                    109: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
                    110: *
                    111: *     .. Scalar Arguments ..
                    112:       INTEGER            INCX, N
                    113:       COMPLEX*16         ALPHA, TAU
                    114: *     ..
                    115: *     .. Array Arguments ..
                    116:       COMPLEX*16         X( * )
                    117: *     ..
                    118: *
                    119: *  =====================================================================
                    120: *
                    121: *     .. Parameters ..
                    122:       DOUBLE PRECISION   ONE, ZERO
                    123:       PARAMETER          ( ONE = 1.0D+0, ZERO = 0.0D+0 )
                    124: *     ..
                    125: *     .. Local Scalars ..
                    126:       INTEGER            J, KNT
                    127:       DOUBLE PRECISION   ALPHI, ALPHR, BETA, RSAFMN, SAFMIN, XNORM
                    128: *     ..
                    129: *     .. External Functions ..
                    130:       DOUBLE PRECISION   DLAMCH, DLAPY3, DZNRM2
                    131:       COMPLEX*16         ZLADIV
                    132:       EXTERNAL           DLAMCH, DLAPY3, DZNRM2, ZLADIV
                    133: *     ..
                    134: *     .. Intrinsic Functions ..
                    135:       INTRINSIC          ABS, DBLE, DCMPLX, DIMAG, SIGN
                    136: *     ..
                    137: *     .. External Subroutines ..
                    138:       EXTERNAL           ZDSCAL, ZSCAL
                    139: *     ..
                    140: *     .. Executable Statements ..
                    141: *
                    142:       IF( N.LE.0 ) THEN
                    143:          TAU = ZERO
                    144:          RETURN
                    145:       END IF
                    146: *
                    147:       XNORM = DZNRM2( N-1, X, INCX )
                    148:       ALPHR = DBLE( ALPHA )
                    149:       ALPHI = DIMAG( ALPHA )
                    150: *
                    151:       IF( XNORM.EQ.ZERO .AND. ALPHI.EQ.ZERO ) THEN
                    152: *
                    153: *        H  =  I
                    154: *
                    155:          TAU = ZERO
                    156:       ELSE
                    157: *
                    158: *        general case
                    159: *
                    160:          BETA = -SIGN( DLAPY3( ALPHR, ALPHI, XNORM ), ALPHR )
                    161:          SAFMIN = DLAMCH( 'S' ) / DLAMCH( 'E' )
                    162:          RSAFMN = ONE / SAFMIN
                    163: *
                    164:          KNT = 0
                    165:          IF( ABS( BETA ).LT.SAFMIN ) THEN
                    166: *
                    167: *           XNORM, BETA may be inaccurate; scale X and recompute them
                    168: *
                    169:    10       CONTINUE
                    170:             KNT = KNT + 1
                    171:             CALL ZDSCAL( N-1, RSAFMN, X, INCX )
                    172:             BETA = BETA*RSAFMN
                    173:             ALPHI = ALPHI*RSAFMN
                    174:             ALPHR = ALPHR*RSAFMN
1.18      bertrand  175:             IF( (ABS( BETA ).LT.SAFMIN) .AND. (KNT .LT. 20) )
1.1       bertrand  176:      $         GO TO 10
                    177: *
                    178: *           New BETA is at most 1, at least SAFMIN
                    179: *
                    180:             XNORM = DZNRM2( N-1, X, INCX )
                    181:             ALPHA = DCMPLX( ALPHR, ALPHI )
                    182:             BETA = -SIGN( DLAPY3( ALPHR, ALPHI, XNORM ), ALPHR )
                    183:          END IF
                    184:          TAU = DCMPLX( ( BETA-ALPHR ) / BETA, -ALPHI / BETA )
                    185:          ALPHA = ZLADIV( DCMPLX( ONE ), ALPHA-BETA )
                    186:          CALL ZSCAL( N-1, ALPHA, X, INCX )
                    187: *
                    188: *        If ALPHA is subnormal, it may lose relative accuracy
                    189: *
                    190:          DO 20 J = 1, KNT
                    191:             BETA = BETA*SAFMIN
                    192:  20      CONTINUE
                    193:          ALPHA = BETA
                    194:       END IF
                    195: *
                    196:       RETURN
                    197: *
                    198: *     End of ZLARFG
                    199: *
                    200:       END

CVSweb interface <joel.bertrand@systella.fr>