LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zhpgv()

subroutine zhpgv ( integer  itype,
character  jobz,
character  uplo,
integer  n,
complex*16, dimension( * )  ap,
complex*16, dimension( * )  bp,
double precision, dimension( * )  w,
complex*16, dimension( ldz, * )  z,
integer  ldz,
complex*16, dimension( * )  work,
double precision, dimension( * )  rwork,
integer  info 
)

ZHPGV

Download ZHPGV + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 ZHPGV computes all the eigenvalues and, optionally, the eigenvectors
 of a complex generalized Hermitian-definite eigenproblem, of the form
 A*x=(lambda)*B*x,  A*Bx=(lambda)*x,  or B*A*x=(lambda)*x.
 Here A and B are assumed to be Hermitian, stored in packed format,
 and B is also positive definite.
Parameters
[in]ITYPE
          ITYPE is INTEGER
          Specifies the problem type to be solved:
          = 1:  A*x = (lambda)*B*x
          = 2:  A*B*x = (lambda)*x
          = 3:  B*A*x = (lambda)*x
[in]JOBZ
          JOBZ is CHARACTER*1
          = 'N':  Compute eigenvalues only;
          = 'V':  Compute eigenvalues and eigenvectors.
[in]UPLO
          UPLO is CHARACTER*1
          = 'U':  Upper triangles of A and B are stored;
          = 'L':  Lower triangles of A and B are stored.
[in]N
          N is INTEGER
          The order of the matrices A and B.  N >= 0.
[in,out]AP
          AP is COMPLEX*16 array, dimension (N*(N+1)/2)
          On entry, the upper or lower triangle of the Hermitian matrix
          A, packed columnwise in a linear array.  The j-th column of A
          is stored in the array AP as follows:
          if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
          if UPLO = 'L', AP(i + (j-1)*(2*n-j)/2) = A(i,j) for j<=i<=n.

          On exit, the contents of AP are destroyed.
[in,out]BP
          BP is COMPLEX*16 array, dimension (N*(N+1)/2)
          On entry, the upper or lower triangle of the Hermitian matrix
          B, packed columnwise in a linear array.  The j-th column of B
          is stored in the array BP as follows:
          if UPLO = 'U', BP(i + (j-1)*j/2) = B(i,j) for 1<=i<=j;
          if UPLO = 'L', BP(i + (j-1)*(2*n-j)/2) = B(i,j) for j<=i<=n.

          On exit, the triangular factor U or L from the Cholesky
          factorization B = U**H*U or B = L*L**H, in the same storage
          format as B.
[out]W
          W is DOUBLE PRECISION array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]Z
          Z is COMPLEX*16 array, dimension (LDZ, N)
          If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
          eigenvectors.  The eigenvectors are normalized as follows:
          if ITYPE = 1 or 2, Z**H*B*Z = I;
          if ITYPE = 3, Z**H*inv(B)*Z = I.
          If JOBZ = 'N', then Z is not referenced.
[in]LDZ
          LDZ is INTEGER
          The leading dimension of the array Z.  LDZ >= 1, and if
          JOBZ = 'V', LDZ >= max(1,N).
[out]WORK
          WORK is COMPLEX*16 array, dimension (max(1, 2*N-1))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (max(1, 3*N-2))
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  ZPPTRF or ZHPEV returned an error code:
             <= N:  if INFO = i, ZHPEV failed to converge;
                    i off-diagonal elements of an intermediate
                    tridiagonal form did not convergeto zero;
             > N:   if INFO = N + i, for 1 <= i <= n, then the leading
                    principal minor of order i of B is not positive.
                    The factorization of B could not be completed and
                    no eigenvalues or eigenvectors were computed.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 163 of file zhpgv.f.

165*
166* -- LAPACK driver routine --
167* -- LAPACK is a software package provided by Univ. of Tennessee, --
168* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
169*
170* .. Scalar Arguments ..
171 CHARACTER JOBZ, UPLO
172 INTEGER INFO, ITYPE, LDZ, N
173* ..
174* .. Array Arguments ..
175 DOUBLE PRECISION RWORK( * ), W( * )
176 COMPLEX*16 AP( * ), BP( * ), WORK( * ), Z( LDZ, * )
177* ..
178*
179* =====================================================================
180*
181* .. Local Scalars ..
182 LOGICAL UPPER, WANTZ
183 CHARACTER TRANS
184 INTEGER J, NEIG
185* ..
186* .. External Functions ..
187 LOGICAL LSAME
188 EXTERNAL lsame
189* ..
190* .. External Subroutines ..
191 EXTERNAL xerbla, zhpev, zhpgst, zpptrf, ztpmv, ztpsv
192* ..
193* .. Executable Statements ..
194*
195* Test the input parameters.
196*
197 wantz = lsame( jobz, 'V' )
198 upper = lsame( uplo, 'U' )
199*
200 info = 0
201 IF( itype.LT.1 .OR. itype.GT.3 ) THEN
202 info = -1
203 ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
204 info = -2
205 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
206 info = -3
207 ELSE IF( n.LT.0 ) THEN
208 info = -4
209 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
210 info = -9
211 END IF
212 IF( info.NE.0 ) THEN
213 CALL xerbla( 'ZHPGV ', -info )
214 RETURN
215 END IF
216*
217* Quick return if possible
218*
219 IF( n.EQ.0 )
220 $ RETURN
221*
222* Form a Cholesky factorization of B.
223*
224 CALL zpptrf( uplo, n, bp, info )
225 IF( info.NE.0 ) THEN
226 info = n + info
227 RETURN
228 END IF
229*
230* Transform problem to standard eigenvalue problem and solve.
231*
232 CALL zhpgst( itype, uplo, n, ap, bp, info )
233 CALL zhpev( jobz, uplo, n, ap, w, z, ldz, work, rwork, info )
234*
235 IF( wantz ) THEN
236*
237* Backtransform eigenvectors to the original problem.
238*
239 neig = n
240 IF( info.GT.0 )
241 $ neig = info - 1
242 IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
243*
244* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
245* backtransform eigenvectors: x = inv(L)**H *y or inv(U)*y
246*
247 IF( upper ) THEN
248 trans = 'N'
249 ELSE
250 trans = 'C'
251 END IF
252*
253 DO 10 j = 1, neig
254 CALL ztpsv( uplo, trans, 'Non-unit', n, bp, z( 1, j ),
255 $ 1 )
256 10 CONTINUE
257*
258 ELSE IF( itype.EQ.3 ) THEN
259*
260* For B*A*x=(lambda)*x;
261* backtransform eigenvectors: x = L*y or U**H *y
262*
263 IF( upper ) THEN
264 trans = 'C'
265 ELSE
266 trans = 'N'
267 END IF
268*
269 DO 20 j = 1, neig
270 CALL ztpmv( uplo, trans, 'Non-unit', n, bp, z( 1, j ),
271 $ 1 )
272 20 CONTINUE
273 END IF
274 END IF
275 RETURN
276*
277* End of ZHPGV
278*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zhpev(jobz, uplo, n, ap, w, z, ldz, work, rwork, info)
ZHPEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices
Definition zhpev.f:138
subroutine zhpgst(itype, uplo, n, ap, bp, info)
ZHPGST
Definition zhpgst.f:113
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine zpptrf(uplo, n, ap, info)
ZPPTRF
Definition zpptrf.f:119
subroutine ztpmv(uplo, trans, diag, n, ap, x, incx)
ZTPMV
Definition ztpmv.f:142
subroutine ztpsv(uplo, trans, diag, n, ap, x, incx)
ZTPSV
Definition ztpsv.f:144
Here is the call graph for this function:
Here is the caller graph for this function: