Go to main content
Oracle Developer Studio 12.5 Man Pages

Exit Print View

Updated: June 2017
 
 

zlaed7 (3p)

Name

zlaed7 - nal matrix after modification by a rank-one symmetric matrix. Used when the original matrix is dense

Synopsis

SUBROUTINE  ZLAED7(  N, CUTPNT, QSIZ, TLVLS, CURLVL, CURPBM, D, Q, LDQ,
RHO, INDXQ,  QSTORE,  QPTR,  PRMPTR,  PERM,  GIVPTR,  GIVCOL,
GIVNUM, WORK, RWORK, IWORK, INFO )


INTEGER CURLVL, CURPBM, CUTPNT, INFO, LDQ, N, QSIZ, TLVLS

DOUBLE PRECISION RHO

INTEGER GIVCOL(2,*), GIVPTR(*), INDXQ(*), IWORK(*), PERM(*), PRMPTR(*),
QPTR(*)

DOUBLE PRECISION D(*), GIVNUM(2,*), QSTORE(*), RWORK(*)

DOUBLE COMPLEX Q(LDQ,*), WORK(*)


SUBROUTINE ZLAED7_64( N, CUTPNT, QSIZ, TLVLS,  CURLVL,  CURPBM,  D,  Q,
LDQ,  RHO, INDXQ, QSTORE, QPTR, PRMPTR, PERM, GIVPTR, GIVCOL,
GIVNUM, WORK, RWORK, IWORK, INFO )


INTEGER*8 CURLVL, CURPBM, CUTPNT, INFO, LDQ, N, QSIZ, TLVLS

DOUBLE PRECISION RHO

INTEGER*8  GIVCOL(2,*),   GIVPTR(*),   INDXQ(*),   IWORK(*),   PERM(*),
PRMPTR(*), QPTR(*)

DOUBLE PRECISION D(*), GIVNUM(2,*), QSTORE(*), RWORK(*)

DOUBLE COMPLEX Q(LDQ,*), WORK(*)


F95 INTERFACE
SUBROUTINE  LAED7(  N,  CUTPNT, QSIZ, TLVLS, CURLVL, CURPBM, D, Q, LDQ,
RHO, INDXQ,  QSTORE,  QPTR,  PRMPTR,  PERM,  GIVPTR,  GIVCOL,
GIVNUM, WORK, RWORK, IWORK, INFO )


INTEGER :: N, CUTPNT, QSIZ, TLVLS, CURLVL, CURPBM, LDQ, INFO

INTEGER, DIMENSION(:) :: INDXQ, QPTR, PRMPTR, PERM, GIVPTR, IWORK

REAL(8), DIMENSION(:,:) :: GIVNUM

COMPLEX(8), DIMENSION(:,:) :: Q

REAL(8), DIMENSION(:) :: D, QSTORE, RWORK

INTEGER, DIMENSION(:,:) :: GIVCOL

COMPLEX(8), DIMENSION(:) :: WORK

REAL(8) :: RHO


SUBROUTINE LAED7_64( N, CUTPNT, QSIZ, TLVLS, CURLVL, CURPBM, D, Q, LDQ,
RHO, INDXQ,  QSTORE,  QPTR,  PRMPTR,  PERM,  GIVPTR,  GIVCOL,
GIVNUM, WORK, RWORK, IWORK, INFO )


INTEGER(8) :: N, CUTPNT, QSIZ, TLVLS, CURLVL, CURPBM, LDQ, INFO

INTEGER(8), DIMENSION(:) :: INDXQ, QPTR, PRMPTR, PERM, GIVPTR, IWORK

REAL(8), DIMENSION(:,:) :: GIVNUM

COMPLEX(8), DIMENSION(:,:) :: Q

REAL(8), DIMENSION(:) :: D, QSTORE, RWORK

INTEGER(8), DIMENSION(:,:) :: GIVCOL

COMPLEX(8), DIMENSION(:) :: WORK

REAL(8) :: RHO


C INTERFACE
#include <sunperf.h>

void  zlaed7  (int  n, int cutpnt, int qsiz, int tlvls, int curlvl, int
curpbm, double *d, doublecomplex *q, int ldq, double rho, int
*indxq,  double  *qstore,  int *qptr, int *prmptr, int *perm,
int *givptr, int *givcol, double *givnum, int *info);


void zlaed7_64 (long n,  long  cutpnt,  long  qsiz,  long  tlvls,  long
curlvl,  long  curpbm, double *d, doublecomplex *q, long ldq,
double rho, long *indxq, double  *qstore,  long  *qptr,  long
*prmptr,  long  *perm,  long  *givptr,  long  *givcol, double
*givnum, long *info);

Description

Oracle Solaris Studio Performance Library                           zlaed7(3P)



NAME
       zlaed7 - is used by sstedc. Compute the updated eigensystem of a diago-
       nal matrix after modification by a rank-one symmetric matrix. Used when
       the original matrix is dense


SYNOPSIS
       SUBROUTINE  ZLAED7(  N, CUTPNT, QSIZ, TLVLS, CURLVL, CURPBM, D, Q, LDQ,
                 RHO, INDXQ,  QSTORE,  QPTR,  PRMPTR,  PERM,  GIVPTR,  GIVCOL,
                 GIVNUM, WORK, RWORK, IWORK, INFO )


       INTEGER CURLVL, CURPBM, CUTPNT, INFO, LDQ, N, QSIZ, TLVLS

       DOUBLE PRECISION RHO

       INTEGER GIVCOL(2,*), GIVPTR(*), INDXQ(*), IWORK(*), PERM(*), PRMPTR(*),
                 QPTR(*)

       DOUBLE PRECISION D(*), GIVNUM(2,*), QSTORE(*), RWORK(*)

       DOUBLE COMPLEX Q(LDQ,*), WORK(*)


       SUBROUTINE ZLAED7_64( N, CUTPNT, QSIZ, TLVLS,  CURLVL,  CURPBM,  D,  Q,
                 LDQ,  RHO, INDXQ, QSTORE, QPTR, PRMPTR, PERM, GIVPTR, GIVCOL,
                 GIVNUM, WORK, RWORK, IWORK, INFO )


       INTEGER*8 CURLVL, CURPBM, CUTPNT, INFO, LDQ, N, QSIZ, TLVLS

       DOUBLE PRECISION RHO

       INTEGER*8  GIVCOL(2,*),   GIVPTR(*),   INDXQ(*),   IWORK(*),   PERM(*),
                 PRMPTR(*), QPTR(*)

       DOUBLE PRECISION D(*), GIVNUM(2,*), QSTORE(*), RWORK(*)

       DOUBLE COMPLEX Q(LDQ,*), WORK(*)


   F95 INTERFACE
       SUBROUTINE  LAED7(  N,  CUTPNT, QSIZ, TLVLS, CURLVL, CURPBM, D, Q, LDQ,
                 RHO, INDXQ,  QSTORE,  QPTR,  PRMPTR,  PERM,  GIVPTR,  GIVCOL,
                 GIVNUM, WORK, RWORK, IWORK, INFO )


       INTEGER :: N, CUTPNT, QSIZ, TLVLS, CURLVL, CURPBM, LDQ, INFO

       INTEGER, DIMENSION(:) :: INDXQ, QPTR, PRMPTR, PERM, GIVPTR, IWORK

       REAL(8), DIMENSION(:,:) :: GIVNUM

       COMPLEX(8), DIMENSION(:,:) :: Q

       REAL(8), DIMENSION(:) :: D, QSTORE, RWORK

       INTEGER, DIMENSION(:,:) :: GIVCOL

       COMPLEX(8), DIMENSION(:) :: WORK

       REAL(8) :: RHO


       SUBROUTINE LAED7_64( N, CUTPNT, QSIZ, TLVLS, CURLVL, CURPBM, D, Q, LDQ,
                 RHO, INDXQ,  QSTORE,  QPTR,  PRMPTR,  PERM,  GIVPTR,  GIVCOL,
                 GIVNUM, WORK, RWORK, IWORK, INFO )


       INTEGER(8) :: N, CUTPNT, QSIZ, TLVLS, CURLVL, CURPBM, LDQ, INFO

       INTEGER(8), DIMENSION(:) :: INDXQ, QPTR, PRMPTR, PERM, GIVPTR, IWORK

       REAL(8), DIMENSION(:,:) :: GIVNUM

       COMPLEX(8), DIMENSION(:,:) :: Q

       REAL(8), DIMENSION(:) :: D, QSTORE, RWORK

       INTEGER(8), DIMENSION(:,:) :: GIVCOL

       COMPLEX(8), DIMENSION(:) :: WORK

       REAL(8) :: RHO


   C INTERFACE
       #include <sunperf.h>

       void  zlaed7  (int  n, int cutpnt, int qsiz, int tlvls, int curlvl, int
                 curpbm, double *d, doublecomplex *q, int ldq, double rho, int
                 *indxq,  double  *qstore,  int *qptr, int *prmptr, int *perm,
                 int *givptr, int *givcol, double *givnum, int *info);


       void zlaed7_64 (long n,  long  cutpnt,  long  qsiz,  long  tlvls,  long
                 curlvl,  long  curpbm, double *d, doublecomplex *q, long ldq,
                 double rho, long *indxq, double  *qstore,  long  *qptr,  long
                 *prmptr,  long  *perm,  long  *givptr,  long  *givcol, double
                 *givnum, long *info);


PURPOSE
       zlaed7 computes the updated eigensystem of a diagonal matrix after mod-
       ification by a rank-one symmetric matrix. This routine is used only for
       the eigenproblem which requires all eigenvalues and  optionally  eigen-
       vectors  of a dense or banded Hermitian matrix that has been reduced to
       tridiagonal form.

       T = Q(in) ( D(in) + RHO *  Z*Z**H  )  Q**H(in)  =  Q(out)  *  D(out)  *
       Q**H(out)

       where  Z = Q**Hu, u is a vector of length N with ones in the CUTPNT and
       CUTPNT + 1 th elements and zeros elsewhere.

       The eigenvectors of the original matrix are stored in Q, and the eigen-
       values are in D.  The algorithm consists of three stages:

       The  first  stage  consists  of  deflating the size of the problem when
       there are multiple eigenvalues or if there is a zero in the  Z  vector.
       For  each  such occurence the dimension of the secular equation problem
       is reduced by one.  This stage is performed by the routine DLAED2.

       The second stage consists of calculating the updated eigenvalues.  This
       is  done  by  finding the roots of the secular equation via the routine
       DLAED4 (as called by SLAED3).  This routine also calculates the  eigen-
       vectors of the current problem.

       The final stage consists of computing the updated eigenvectors directly
       using the updated eigenvalues.  The eigenvectors for the current  prob-
       lem are multiplied with the eigenvectors from the overall problem.


ARGUMENTS
       N (input)
                 N is INTEGER
                 The dimension of the symmetric tridiagonal matrix.  N >= 0.


       CUTPNT (input)
                 CUTPNT is INTEGER
                 Contains the location of the last eigenvalue in the leading
                 sub-matrix.  min(1,N) <= CUTPNT <= N.


       QSIZ (input)
                 QSIZ is INTEGER
                 The dimension of the unitary matrix used to reduce
                 the full matrix to tridiagonal form.  QSIZ >= N.


       TLVLS (input)
                 TLVLS is INTEGER
                 The total number of merging levels in the overall divide and
                 conquer tree.


       CURLVL (input)
                 CURLVL is INTEGER
                 The current level in the overall merge routine,
                 0 <= curlvl <= tlvls.


       CURPBM (input)
                 CURPBM is INTEGER
                 The current problem in the current level in the overall
                 merge routine (counting from upper left to lower right).


       D (input/output)
                 D is DOUBLE PRECISION array, dimension (N)
                 On entry, the eigenvalues of the rank-1-perturbed matrix.
                 On exit, the eigenvalues of the repaired matrix.


       Q (input/output)
                 Q is COMPLEX*16 array, dimension (LDQ,N)
                 On entry, the eigenvectors of the rank-1-perturbed matrix.
                 On exit, the eigenvectors of the repaired tridiagonal matrix.


       LDQ (input)
                 LDQ is INTEGER
                 The leading dimension of the array Q.  LDQ >= max(1,N).


       RHO (input)
                 RHO is DOUBLE PRECISION
                 Contains the subdiagonal element used to create the rank-1
                 modification.


       INDXQ (output)
                 INDXQ is INTEGER array, dimension (N)
                 This contains the permutation which will reintegrate the
                 subproblem just solved back into sorted order,
                 ie. D( INDXQ( I = 1, N ) ) will be in ascending order.


       IWORK (output)
                 IWORK is INTEGER array, dimension (4*N)


       RWORK (output)
                 RWORK is DOUBLE PRECISION array,
                 dimension (3*N+2*QSIZ*N)


       WORK (output)
                 WORK is COMPLEX*16 array, dimension (QSIZ*N)


       QSTORE (input/output)
                 QSTORE is DOUBLE PRECISION array, dimension (N**2+1)
                 Stores eigenvectors of submatrices encountered during
                 divide and conquer, packed together. QPTR points to
                 beginning of the submatrices.


       QPTR (input/output)
                 QPTR is INTEGER array, dimension (N+2)
                 List of indices pointing to beginning of submatrices stored
                 in QSTORE. The submatrices are numbered starting at the
                 bottom left of the divide and conquer tree, from left to
                 right and bottom to top.


       PRMPTR (input)
                 PRMPTR is INTEGER array, dimension (N lg N)
                 Contains a list of pointers which indicate where in PERM a
                 level's permutation is stored.  PRMPTR(i+1) - PRMPTR(i)
                 indicates the size of the permutation and also the size of
                 the full, non-deflated problem.


       PERM (input)
                 PERM is INTEGER array, dimension (N lg N)
                 Contains the permutations (from deflation and sorting) to be
                 applied to each eigenblock.


       GIVPTR (input)
                 GIVPTR is INTEGER array, dimension (N lg N)
                 Contains a list of pointers which indicate where in GIVCOL a
                 level's Givens rotations are stored.  GIVPTR(i+1) - GIVPTR(i)
                 indicates the number of Givens rotations.


       GIVCOL (input)
                 GIVCOL is INTEGER array, dimension (2, N lg N)
                 Each  pair  of  numbers  indicates  a pair of columns to take
                 place
                 in a Givens rotation.


       GIVNUM (input)
                 GIVNUM is DOUBLE PRECISION array, dimension (2, N lg N)
                 Each number indicates the S value to be used in the
                 corresponding Givens rotation.


       INFO (output)
                 INFO is INTEGER
                 = 0:  successful exit.
                 < 0:  if INFO = -i, the i-th argument had an illegal value.
                 > 0:  if INFO = 1, an eigenvalue did not converge




                                  7 Nov 2015                        zlaed7(3P)