PDLAQR2(3)    ScaLAPACK routine of NEC Numeric Library Collection   PDLAQR2(3)



NAME
       PDLAQR2 - accepts as input an upper Hessenberg matrix A and performs an
       orthogonal similarity transformation designed  to  detect  and  deflate
       fully converged eigenvalues from a trailing principal submatrix

SYNOPSIS
       SUBROUTINE PDLAQR2( WANTT,  WANTZ,  N,  KTOP, KBOT, NW, A, DESCA, ILOZ,
                           IHIZ, Z, DESCZ, NS, ND, SR, SI, T, LDT, V, LDV, WR,
                           WI, WORK, LWORK )

           INTEGER         IHIZ, ILOZ, KBOT, KTOP, LDT, LDV, LWORK, N, ND, NS,
                           NW

           LOGICAL         WANTT, WANTZ

           INTEGER         DESCA( * ), DESCZ( * )

           DOUBLE          PRECISION A( * ), SI( KBOT ), SR( KBOT ), T( LDT, *
                           ), V( LDV, * ), WORK( * ), WI( * ), WR( * ), Z( * )

PURPOSE
       Aggressive early deflation:

       PDLAQR2 accepts as input an upper Hessenberg matrix A and  performs  an
       orthogonal  similarity  transformation  designed  to detect and deflate
       fully converged eigenvalues from a trailing  principal  submatrix.   On
       output A has been overwritten by a new Hessenberg matrix that is a per-
       turbation of an orthogonal similarity transformation of A.  It is to be
       hoped that the final version of H has many zero subdiagonal entries.

       This routine handles small deflation windows which is affordable by one
       processor. Normally, it is  called  by  PDLAQR1.  All  the  inputs  are
       assumed to be valid without checking.


       Notes
       =====

       Each  global data object is described by an associated description vec-
       tor.  This vector stores the information required to establish the map-
       ping between an object element and its corresponding process and memory
       location.

       Let A be a generic term for any 2D block  cyclicly  distributed  array.
       Such a global array has an associated description vector DESCA.  In the
       following comments, the character _ should be read as  "of  the  global
       array".

       NOTATION        STORED IN      EXPLANATION
       --------------- -------------- --------------------------------------
       DTYPE_A(global) DESCA( DTYPE_ )The descriptor type.  In this case,
                                      DTYPE_A = 1.
       CTXT_A (global) DESCA( CTXT_ ) The BLACS context handle, indicating
                                      the BLACS process grid A is distribu-
                                      ted over. The context itself is glo-
                                      bal, but the handle (the integer
                                      value) may vary.
       M_A    (global) DESCA( M_ )    The number of rows in the global
                                      array A.
       N_A    (global) DESCA( N_ )    The number of columns in the global
                                      array A.
       MB_A   (global) DESCA( MB_ )   The blocking factor used to distribute
                                      the rows of the array.
       NB_A   (global) DESCA( NB_ )   The blocking factor used to distribute
                                      the columns of the array.
       RSRC_A (global) DESCA( RSRC_ ) The process row over which the first
                                      row of the array A is distributed.
       CSRC_A (global) DESCA( CSRC_ ) The process column over which the
                                      first column of the array A is
                                      distributed.
       LLD_A  (local)  DESCA( LLD_ )  The leading dimension of the local
                                      array.  LLD_A >= MAX(1,LOCr(M_A)).

       Let  K  be  the  number of rows or columns of a distributed matrix, and
       assume that its process grid has dimension p x q.
       LOCr( K ) denotes the number of elements of  K  that  a  process  would
       receive  if K were distributed over the p processes of its process col-
       umn.
       Similarly, LOCc( K ) denotes the number of elements of K that a process
       would receive if K were distributed over the q processes of its process
       row.
       The values of LOCr() and LOCc() may be determined via  a  call  to  the
       ScaLAPACK tool function, NUMROC:
               LOCr( M ) = NUMROC( M, MB_A, MYROW, RSRC_A, NPROW ),
               LOCc( N ) = NUMROC( N, NB_A, MYCOL, CSRC_A, NPCOL ).
       An upper bound for these quantities may be computed by:
               LOCr( M ) <= ceil( ceil(M/MB_A)/NPROW )*MB_A
               LOCc( N ) <= ceil( ceil(N/NB_A)/NPCOL )*NB_A


ARGUMENTS
       WANTT   (global input) LOGICAL
               If  .TRUE.,  then  the  Hessenberg matrix H is fully updated so
               that the quasi-triangular Schur  factor  may  be  computed  (in
               cooperation with the calling subroutine).
               If  .FALSE.,  then  only enough of H is updated to preserve the
               eigenvalues.

       WANTZ   (global input) LOGICAL
               If .TRUE., then the orthogonal matrix Z is updated so  so  that
               the  orthogonal  Schur  factor  may be computed (in cooperation
               with the calling subroutine).
               If .FALSE., then Z is not referenced.

       N       (global input) INTEGER
               The order of the matrix H and (if WANTZ is .TRUE.) the order of
               the orthogonal matrix Z.

       KTOP    (global input) INTEGER

       KBOT    (global input) INTEGER
               It  is  assumed  without  a  check  that  either  KBOT  =  N or
               H(KBOT+1,KBOT)=0.  KBOT and KTOP together determine an isolated
               block  along  the  diagonal  of the Hessenberg matrix. However,
               H(KTOP,KTOP-1)=0 is  not  essentially  necessary  if  WANTT  is
               .TRUE. .

       NW      (global input) INTEGER
               Deflation window size.  1 .LE. NW .LE. (KBOT-KTOP+1).  Normally
               NW .GE. 3 if PDLAQR2 is called by PDLAQR1.

       A       (local input/output) DOUBLE PRECISION array, dimension
               (DESCH(LLD_),*) On input the initial N-by-N section of A stores
               the Hessenberg matrix undergoing aggressive early deflation.
               On  output  A  has been transformed by an orthogonal similarity
               transformation, perturbed, and the returned to Hessenberg  form
               that (it is to be hoped) has some zero subdiagonal entries.

       DESCA   (global and local input) INTEGER array of dimension DLEN_.
               The array descriptor for the distributed matrix A.

       ILOZ    (global input) INTEGER

       IHIZ    (global input) INTEGER
               Specify  the rows of Z to which transformations must be applied
               if WANTZ is .TRUE.. 1 .LE. ILOZ .LE. IHIZ .LE. N.

       Z       (input/output) DOUBLE PRECISION array, dimension
               (DESCH(LLD_),*)
               IF WANTZ is .TRUE., then on output, the  orthogonal  similarity
               transformation   mentioned  above  has  been  accumulated  into
               Z(ILOZ:IHIZ,ILO:IHI) from the right.
               If WANTZ is .FALSE., then Z is unreferenced.

       DESCZ   (global and local input) INTEGER array of dimension DLEN_.
               The array descriptor for the distributed matrix Z.

       NS      (global output) INTEGER
               The number of unconverged (ie approximate) eigenvalues returned
               in  SR and SI that may be used as shifts by the calling subrou-
               tine.

       ND      (global output) INTEGER
               The number of converged eigenvalues uncovered by  this  subrou-
               tine.

       SR      (global output) DOUBLE PRECISION array, dimension KBOT

       SI      (global output) DOUBLE PRECISION array, dimension KBOT
               On  output,  the real and imaginary parts of approximate eigen-
               values that may be used for shifts are  stored  in  SR(KBOT-ND-
               NS+1) through SR(KBOT-ND) and SI(KBOT-ND-NS+1) through SI(KBOT-
               ND), respectively.
               On proc #0, the real and imaginary parts of converged eigenval-
               ues  are  stored in SR(KBOT-ND+1) through SR(KBOT) and SI(KBOT-
               ND+1) through  SI(KBOT),  respectively.  On  other  processors,
               these entries are set to zero.

       T       (local workspace) DOUBLE PRECISION array, dimension LDT*NW.

       LDT     (local input) INTEGER
               The leading dimension of the array T.
               LDT >= NW.

       V       (local workspace) DOUBLE PRECISION array, dimension LDV*NW.

       LDV     (local input) INTEGER
               The leading dimension of the array V.
               LDV >= NW.

       WR      (local workspace) DOUBLE PRECISION array, dimension KBOT.

       WI      (local workspace) DOUBLE PRECISION array, dimension KBOT.

       WORK    (local workspace) DOUBLE PRECISION array, dimension LWORK.

       LWORK   (local input) INTEGER
               WORK(LWORK) is a local array and LWORK is assumed big enough so
               that LWORK >= NW*NW.



ScaLAPACK routine               31 October 2017                     PDLAQR2(3)