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

◆ chpgvd()

subroutine chpgvd ( integer  itype,
character  jobz,
character  uplo,
integer  n,
complex, dimension( * )  ap,
complex, dimension( * )  bp,
real, dimension( * )  w,
complex, dimension( ldz, * )  z,
integer  ldz,
complex, dimension( * )  work,
integer  lwork,
real, dimension( * )  rwork,
integer  lrwork,
integer, dimension( * )  iwork,
integer  liwork,
integer  info 
)

CHPGVD

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

Purpose:
 CHPGVD 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.
 If eigenvectors are desired, it uses a divide and conquer algorithm.
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 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 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 REAL array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]Z
          Z is COMPLEX 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 array, dimension (MAX(1,LWORK))
          On exit, if INFO = 0, WORK(1) returns the required LWORK.
[in]LWORK
          LWORK is INTEGER
          The dimension of array WORK.
          If N <= 1,               LWORK >= 1.
          If JOBZ = 'N' and N > 1, LWORK >= N.
          If JOBZ = 'V' and N > 1, LWORK >= 2*N.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the required sizes of the WORK, RWORK and
          IWORK arrays, returns these values as the first entries of
          the WORK, RWORK and IWORK arrays, and no error message
          related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]RWORK
          RWORK is REAL array, dimension (MAX(1,LRWORK))
          On exit, if INFO = 0, RWORK(1) returns the required LRWORK.
[in]LRWORK
          LRWORK is INTEGER
          The dimension of array RWORK.
          If N <= 1,               LRWORK >= 1.
          If JOBZ = 'N' and N > 1, LRWORK >= N.
          If JOBZ = 'V' and N > 1, LRWORK >= 1 + 5*N + 2*N**2.

          If LRWORK = -1, then a workspace query is assumed; the
          routine only calculates the required sizes of the WORK, RWORK
          and IWORK arrays, returns these values as the first entries
          of the WORK, RWORK and IWORK arrays, and no error message
          related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]IWORK
          IWORK is INTEGER array, dimension (MAX(1,LIWORK))
          On exit, if INFO = 0, IWORK(1) returns the required LIWORK.
[in]LIWORK
          LIWORK is INTEGER
          The dimension of array IWORK.
          If JOBZ  = 'N' or N <= 1, LIWORK >= 1.
          If JOBZ  = 'V' and N > 1, LIWORK >= 3 + 5*N.

          If LIWORK = -1, then a workspace query is assumed; the
          routine only calculates the required sizes of the WORK, RWORK
          and IWORK arrays, returns these values as the first entries
          of the WORK, RWORK and IWORK arrays, and no error message
          related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  CPPTRF or CHPEVD returned an error code:
             <= N:  if INFO = i, CHPEVD 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.
Contributors:
Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA

Definition at line 223 of file chpgvd.f.

225*
226* -- LAPACK driver routine --
227* -- LAPACK is a software package provided by Univ. of Tennessee, --
228* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
229*
230* .. Scalar Arguments ..
231 CHARACTER JOBZ, UPLO
232 INTEGER INFO, ITYPE, LDZ, LIWORK, LRWORK, LWORK, N
233* ..
234* .. Array Arguments ..
235 INTEGER IWORK( * )
236 REAL RWORK( * ), W( * )
237 COMPLEX AP( * ), BP( * ), WORK( * ), Z( LDZ, * )
238* ..
239*
240* =====================================================================
241*
242* .. Local Scalars ..
243 LOGICAL LQUERY, UPPER, WANTZ
244 CHARACTER TRANS
245 INTEGER J, LIWMIN, LRWMIN, LWMIN, NEIG
246* ..
247* .. External Functions ..
248 LOGICAL LSAME
249 REAL SROUNDUP_LWORK
250 EXTERNAL lsame, sroundup_lwork
251* ..
252* .. External Subroutines ..
253 EXTERNAL chpevd, chpgst, cpptrf, ctpmv, ctpsv, xerbla
254* ..
255* .. Intrinsic Functions ..
256 INTRINSIC max, real
257* ..
258* .. Executable Statements ..
259*
260* Test the input parameters.
261*
262 wantz = lsame( jobz, 'V' )
263 upper = lsame( uplo, 'U' )
264 lquery = ( lwork.EQ.-1 .OR. lrwork.EQ.-1 .OR. liwork.EQ.-1 )
265*
266 info = 0
267 IF( itype.LT.1 .OR. itype.GT.3 ) THEN
268 info = -1
269 ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
270 info = -2
271 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
272 info = -3
273 ELSE IF( n.LT.0 ) THEN
274 info = -4
275 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
276 info = -9
277 END IF
278*
279 IF( info.EQ.0 ) THEN
280 IF( n.LE.1 ) THEN
281 lwmin = 1
282 liwmin = 1
283 lrwmin = 1
284 ELSE
285 IF( wantz ) THEN
286 lwmin = 2*n
287 lrwmin = 1 + 5*n + 2*n**2
288 liwmin = 3 + 5*n
289 ELSE
290 lwmin = n
291 lrwmin = n
292 liwmin = 1
293 END IF
294 END IF
295*
296 work( 1 ) = sroundup_lwork(lwmin)
297 rwork( 1 ) = lrwmin
298 iwork( 1 ) = liwmin
299 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
300 info = -11
301 ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
302 info = -13
303 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
304 info = -15
305 END IF
306 END IF
307*
308 IF( info.NE.0 ) THEN
309 CALL xerbla( 'CHPGVD', -info )
310 RETURN
311 ELSE IF( lquery ) THEN
312 RETURN
313 END IF
314*
315* Quick return if possible
316*
317 IF( n.EQ.0 )
318 $ RETURN
319*
320* Form a Cholesky factorization of B.
321*
322 CALL cpptrf( uplo, n, bp, info )
323 IF( info.NE.0 ) THEN
324 info = n + info
325 RETURN
326 END IF
327*
328* Transform problem to standard eigenvalue problem and solve.
329*
330 CALL chpgst( itype, uplo, n, ap, bp, info )
331 CALL chpevd( jobz, uplo, n, ap, w, z, ldz, work, lwork, rwork,
332 $ lrwork, iwork, liwork, info )
333 lwmin = int( max( real( lwmin ), real( work( 1 ) ) ) )
334 lrwmin = int( max( real( lrwmin ), real( rwork( 1 ) ) ) )
335 liwmin = int( max( real( liwmin ), real( iwork( 1 ) ) ) )
336*
337 IF( wantz ) THEN
338*
339* Backtransform eigenvectors to the original problem.
340*
341 neig = n
342 IF( info.GT.0 )
343 $ neig = info - 1
344 IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
345*
346* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
347* backtransform eigenvectors: x = inv(L)**H *y or inv(U)*y
348*
349 IF( upper ) THEN
350 trans = 'N'
351 ELSE
352 trans = 'C'
353 END IF
354*
355 DO 10 j = 1, neig
356 CALL ctpsv( uplo, trans, 'Non-unit', n, bp, z( 1, j ),
357 $ 1 )
358 10 CONTINUE
359*
360 ELSE IF( itype.EQ.3 ) THEN
361*
362* For B*A*x=(lambda)*x;
363* backtransform eigenvectors: x = L*y or U**H *y
364*
365 IF( upper ) THEN
366 trans = 'C'
367 ELSE
368 trans = 'N'
369 END IF
370*
371 DO 20 j = 1, neig
372 CALL ctpmv( uplo, trans, 'Non-unit', n, bp, z( 1, j ),
373 $ 1 )
374 20 CONTINUE
375 END IF
376 END IF
377*
378 work( 1 ) = sroundup_lwork(lwmin)
379 rwork( 1 ) = lrwmin
380 iwork( 1 ) = liwmin
381 RETURN
382*
383* End of CHPGVD
384*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine chpevd(jobz, uplo, n, ap, w, z, ldz, work, lwork, rwork, lrwork, iwork, liwork, info)
CHPEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition chpevd.f:194
subroutine chpgst(itype, uplo, n, ap, bp, info)
CHPGST
Definition chpgst.f:113
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine cpptrf(uplo, n, ap, info)
CPPTRF
Definition cpptrf.f:119
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
subroutine ctpmv(uplo, trans, diag, n, ap, x, incx)
CTPMV
Definition ctpmv.f:142
subroutine ctpsv(uplo, trans, diag, n, ap, x, incx)
CTPSV
Definition ctpsv.f:144
Here is the call graph for this function:
Here is the caller graph for this function: