Developer Reference for Intel® oneAPI Math Kernel Library for Fortran

ID 766686
Date 12/16/2022
Public

A newer version of this document is available. Customers should click here to go to the newest version.

Document Table of Contents

p?hegvx

Computes selected eigenvalues and, optionally, eigenvectors of a complex generalized Hermitian positive-definite eigenproblem.

Syntax

call pchegvx(ibtype, jobz, range, uplo, n, a, ia, ja, desca, b, ib, jb, descb, vl, vu, il, iu, abstol, m, nz, w, orfac, z, iz, jz, descz, work, lwork, rwork, lrwork, iwork, liwork, ifail, iclustr, gap, info)

call pzhegvx(ibtype, jobz, range, uplo, n, a, ia, ja, desca, b, ib, jb, descb, vl, vu, il, iu, abstol, m, nz, w, orfac, z, iz, jz, descz, work, lwork, rwork, lrwork, iwork, liwork, ifail, iclustr, gap, info)

Include Files

Description

The p?hegvx routine computes all the eigenvalues, and optionally, the eigenvectors of a complex generalized Hermitian positive-definite eigenproblem, of the form

sub(A)*x = λ*sub(B)*x, sub(A)*sub(B)*x = λ*x, or sub(B)*sub(A)*x = λ*x.

Here sub (A) denoting A(ia:ia+n-1, ja:ja+n-1) and sub(B) are assumed to be Hermitian and sub(B) denoting B(ib:ib+n-1, jb:jb+n-1) is also positive definite.

Product and Performance Information

Performance varies by use, configuration and other factors. Learn more at www.Intel.com/PerformanceIndex.

Notice revision #20201201

Input Parameters
ibtype

(global) INTEGER. Must be 1 or 2 or 3.

Specifies the problem type to be solved:

If ibtype = 1, the problem type is

sub(A)*x = lambda*sub(B)*x;

If ibtype = 2, the problem type is

sub(A)*sub(B)*x = lambda*x;

If ibtype = 3, the problem type is

sub(B)*sub(A)*x = lambda*x.

jobz

(global) CHARACTER*1. Must be 'N' or 'V'.

If jobz ='N', then compute eigenvalues only.

If jobz ='V', then compute eigenvalues and eigenvectors.

range

(global) CHARACTER*1. Must be 'A' or 'V' or 'I'.

If range = 'A', the routine computes all eigenvalues.

If range = 'V', the routine computes eigenvalues in the interval: [vl, vu]

If range = 'I', the routine computes eigenvalues with indices il through iu.

uplo

(global) CHARACTER*1. Must be 'U' or 'L'.

If uplo = 'U', arrays a and b store the upper triangles of sub(A) and sub (B);

If uplo = 'L', arrays a and b store the lower triangles of sub(A) and sub (B).

n

(global) INTEGER.

The order of the matrices sub(A) and sub (B) (n 0).

a

(local)

COMPLEX for pchegvx

DOUBLE COMPLEX for pzhegvx.

Pointer into the local memory to an array of size (lld_a, LOCc(ja+n-1)). On entry, this array contains the local pieces of the n-by-n Hermitian distributed matrix sub(A). If uplo = 'U', the leading n-by-n upper triangular part of sub(A) contains the upper triangular part of the matrix. If uplo = 'L', the leading n-by-n lower triangular part of sub(A) contains the lower triangular part of the matrix.

ia, ja

(global) INTEGER.

The row and column indices in the global matrix A indicating the first row and the first column of the submatrix A, respectively.

desca

(global and local) INTEGER array of size dlen_.

The array descriptor for the distributed matrix A. If desca(ctxt_) is incorrect, p?hegvx cannot guarantee correct error reporting.

b

(local).

COMPLEX for pchegvx

DOUBLE COMPLEX for pzhegvx.

Pointer into the local memory to an array of size (lld_b, LOCc(jb+n-1)). On entry, this array contains the local pieces of the n-by-n Hermitian distributed matrix sub(B).

If uplo = 'U', the leading n-by-n upper triangular part of sub(B) contains the upper triangular part of the matrix.

If uplo = 'L', the leading n-by-n lower triangular part of sub(B) contains the lower triangular part of the matrix.

ib, jb

(global) INTEGER.

The row and column indices in the global matrix B indicating the first row and the first column of the submatrix B, respectively.

descb

(global and local) INTEGER array of size dlen_.

The array descriptor for the distributed matrix B.descb(ctxt_) must be equal to desca(ctxt_).

vl, vu

(global)

REAL for pchegvx

DOUBLE PRECISION for pzhegvx.

If range = 'V', the lower and upper bounds of the interval to be searched for eigenvalues.

If range = 'A' or 'I', vl and vu are not referenced.

il, iu

(global)

INTEGER.

If range = 'I', the indices in ascending order of the smallest and largest eigenvalues to be returned. Constraint: il 1, min(il, n) ≤ iun

If range = 'A' or 'V', il and iu are not referenced.

abstol

(global)

REAL for pchegvx

DOUBLE PRECISION for pzhegvx.

If jobz='V', setting abstol to p?lamch(context, 'U') yields the most orthogonal eigenvectors.

The absolute error tolerance for the eigenvalues. An approximate eigenvalue is accepted as converged when it is determined to lie in an interval [a,b] of width less than or equal to

abstol + eps*max(|a|,|b|),

where eps is the machine precision. If abstol is less than or equal to zero, then eps*norm(T) will be used in its place, where norm(T) is the 1-norm of the tridiagonal matrix obtained by reducing A to tridiagonal form.

Eigenvalues will be computed most accurately when abstol is set to twice the underflow threshold 2*p?lamch('S') not zero. If this routine returns with ((mod(info,2)0).or. * (mod(info/8,2)0)), indicating that some eigenvalues or eigenvectors did not converge, try setting abstol to 2*p?lamch('S').

NOTE:

mod(x,y) is the integer remainder of x/y.

orfac

(global).

REAL for pchegvx

DOUBLE PRECISION for pzhegvx.

Specifies which eigenvectors should be reorthogonalized. Eigenvectors that correspond to eigenvalues which are within tol=orfac*norm(A) of each other are to be reorthogonalized. However, if the workspace is insufficient (see lwork), tol may be decreased until all eigenvectors to be reorthogonalized can be stored in one process. No reorthogonalization will be done if orfac equals zero. A default value of 1.0E-3 is used if orfac is negative. orfac should be identical on all processes.

iz, jz

(global) INTEGER. The row and column indices in the global matrix Z indicating the first row and the first column of the submatrix Z, respectively.

descz

(global and local) INTEGER array of size dlen_. The array descriptor for the distributed matrix Z.descz( ctxt_ ) must equal desca( ctxt_ ).

work

(local)

COMPLEX for pchegvx

DOUBLE COMPLEX for pzhegvx.

Workspace array of size lwork

lwork

(local).

INTEGER. The size of the array work.

If only eigenvalues are requested:

lwork n+ max(NB*(np0 + 1), 3)

If eigenvectors are requested:

lworkn + (np0+ mq0 + NB)*NB

with nq0 = numroc(nn, NB, 0, 0, NPCOL).

For optimal performance, greater workspace is needed, that is

lworkmax(lwork, n, nhetrd_lwopt, nhegst_lwopt)

where lwork is as defined above, and

nhetrd_lwork = 2*(anb+1)*(4*nps+2) + (nps + 1)*nps;

nhegst_lwopt = 2*np0*nb + nq0*nb + nb*nb

nb = desca(mb_)

np0 = numroc(n, nb, 0, 0, NPROW)

nq0 = numroc(n, nb, 0, 0, NPCOL)

ictxt = desca(ctxt_)

anb = pjlaenv(ictxt, 3, 'p?hettrd', 'L', 0, 0, 0, 0)

sqnpc = sqrt(dble(NPROW * NPCOL))

nps = max(numroc(n, 1, 0, 0, sqnpc), 2*anb)

numroc is a ScaLAPACK tool functions;

pjlaenv is a ScaLAPACK environmental inquiry function MYROW, MYCOL, NPROW and NPCOL can be determined by calling the subroutine blacs_gridinfo.

If lwork = -1, then lwork is global input and a workspace query is assumed; the routine only calculates the size required for optimal performance for all work arrays. Each of these values is returned in the first entry of the corresponding work arrays, and no error message is issued by pxerbla.

rwork

(local)

REAL for pchegvx

DOUBLE PRECISION for pzhegvx.

Workspace array of size lrwork.

lrwork

(local) INTEGER. The size of the array rwork.

See below for definitions of variables used to define lrwork.

If no eigenvectors are requested (jobz = 'N'), then lrwork ≥ 5*nn+4*n

If eigenvectors are requested (jobz = 'V'), then the amount of workspace required to guarantee that all eigenvectors are computed is:

lrwork ≥ 4*n + max(5*nn, np0*mq0)+iceil(neig, NPROW*NPCOL)*nn

The computed eigenvectors may not be orthogonal if the minimal workspace is supplied and orfac is too small. If you want to guarantee orthogonality (at the cost of potentially poor performance) you should add the following value to lrwork:

(clustersize-1)*n,

where clustersize is the number of eigenvalues in the largest cluster, where a cluster is defined as a set of close eigenvalues:

{w(k),..., w(k+clustersize-1)|w(j+1) ≤ w(j)+orfac*2*norm(A)}

Variable definitions:

neig = number of eigenvectors requested;

nb = desca(mb_) = desca(nb_) = descz(mb_) = descz(nb_);

nn = max(n, nb, 2);

desca(rsrc_) = desca(nb_) = descz(rsrc_) = descz(csrc_) = 0 ;

np0 = numroc(nn, nb, 0, 0, NPROW);

mq0 = numroc(max(neig, nb, 2), nb, 0, 0, NPCOL);

iceil(x, y) is a ScaLAPACK function returning ceiling(x/y).

When lrwork is too small:

If lwork is too small to guarantee orthogonality, p?hegvx attempts to maintain orthogonality in the clusters with the smallest spacing between the eigenvalues.

If lwork is too small to compute all the eigenvectors requested, no computation is performed and info= -25 is returned. Note that when range='V', p?hegvx does not know how many eigenvectors are requested until the eigenvalues are computed. Therefore, when range='V' and as long as lwork is large enough to allow p?hegvx to compute the eigenvalues, p?hegvx will compute the eigenvalues and as many eigenvectors as it can.

Relationship between workspace, orthogonality & performance:

If clustersize > n/sqrt(NPROW*NPCOL), then providing enough space to compute all the eigenvectors orthogonally will cause serious degradation in performance. In the limit (that is, clustersize = n-1) p?stein will perform no better than ?stein on 1 processor.

For clustersize = n/sqrt(NPROW*NPCOL) reorthogonalizing all eigenvectors will increase the total execution time by a factor of 2 or more.

For clustersize>n/sqrt(NPROW*NPCOL) execution time will grow as the square of the cluster size, all other factors remaining equal and assuming enough workspace. Less workspace means less reorthogonalization but faster execution.

If lwork = -1, then lrwork is global input and a workspace query is assumed; the routine only calculates the size required for optimal performance for all work arrays. Each of these values is returned in the first entry of the corresponding work arrays, and no error message is issued by pxerbla.

iwork

(local) INTEGER. Workspace array.

liwork

(local) INTEGER , size of iwork.

liwork ≥ 6*nnp

Where: nnp = max(n, NPROW*NPCOL + 1, 4)

If liwork = -1, then liwork is global input and a workspace query is assumed; the routine only calculates the minimum and optimal size for all work arrays. Each of these values is returned in the first entry of the corresponding work array, and no error message is issued by pxerbla.

Output Parameters
a

On exit, if jobz = 'V', then if info = 0, sub(A) contains the distributed matrix Z of eigenvectors.

The eigenvectors are normalized as follows:

If ibtype = 1 or 2, then ZH*sub(B)*Z = i;

If ibtype = 3, then ZH*inv(sub(B))*Z = i.

If jobz = 'N', then on exit the upper triangle (if uplo='U') or the lower triangle (if uplo='L') of sub(A), including the diagonal, is destroyed.

b

On exit, if infon, the part of sub(B) containing the matrix is overwritten by the triangular factor U or L from the Cholesky factorization sub(B) = UH*U, or sub(B) = L*LH.

m

(global) INTEGER. The total number of eigenvalues found, 0 ≤ mn.

nz

(global) INTEGER. Total number of eigenvectors computed. 0 < nz < m. The number of columns of z that are filled.

If jobz'V', nz is not referenced.

If jobz = 'V', nz = m unless the user supplies insufficient space and p?hegvx is not able to detect this before beginning computation. To get all the eigenvectors requested, the user must supply both sufficient space to hold the eigenvectors in z (mdescz(n_)) and sufficient workspace to compute them. (See lwork below.) The routine p?hegvx is always able to detect insufficient space without computation unless range = 'V'.

w

(global)

REAL for pchegvx

DOUBLE PRECISION for pzhegvx.

Array of size n. On normal exit, the first m entries contain the selected eigenvalues in ascending order.

z

(local).

COMPLEX for pchegvx

DOUBLE COMPLEX for pzhegvx.

global size (n, n), local size (lld_z, LOCc(jz+n-1)).

If jobz = 'V', then on normal exit the first m columns of z contain the orthonormal eigenvectors of the matrix corresponding to the selected eigenvalues. If an eigenvector fails to converge, then that column of z contains the latest approximation to the eigenvector, and the index of the eigenvector is returned in ifail.

If jobz = 'N', then z is not referenced.

work

On exit, work(1) returns the optimal amount of workspace.

rwork

On exit, rwork(1) contains the amount of workspace required for optimal efficiency

If jobz='N'rwork(1) = optimal amount of workspace required to compute eigenvalues efficiently

If jobz='V'rwork(1) = optimal amount of workspace required to compute eigenvalues and eigenvectors efficiently with no guarantee on orthogonality.

If range='V', it is assumed that all eigenvectors may be required when computing optimal workspace.

ifail

(global) INTEGER.

Array of size n.

ifail provides additional information when info0

If (mod(info/16,2)0), then ifail(1) indicates the order of the smallest minor which is not positive definite.

If (mod(info,2)0) on exit, then ifail(1) contains the indices of the eigenvectors that failed to converge.

If neither of the above error conditions are held, and jobz = 'V', then the first m elements of ifail are set to zero.

iclustr

(global) INTEGER.

Array of size (2*NPROW*NPCOL). This array contains indices of eigenvectors corresponding to a cluster of eigenvalues that could not be reorthogonalized due to insufficient workspace (see lwork, orfac and info). Eigenvectors corresponding to clusters of eigenvalues indexed iclustr(2*i-1) to iclustr(2*i), could not be reorthogonalized due to lack of workspace. Hence the eigenvectors corresponding to these clusters may not be orthogonal.

iclustr() is a zero terminated array. (iclustr(2*k)0.and.clustr(2*k+1)=0) if and only if k is the number of clusters.

iclustr is not referenced if jobz = 'N'.

gap

(global)

REAL for pchegvx

DOUBLE PRECISION for pzhegvx.

Array of size NPROW*NPCOL.

This array contains the gap between eigenvalues whose eigenvectors could not be reorthogonalized. The output values in this array correspond to the clusters indicated by the array iclustr. As a result, the dot product between eigenvectors corresponding to the i-th cluster may be as high as (C*n)/gap(i), where C is a small constant.

info

(global) INTEGER.

If info = 0, the execution is successful.

If info <0: the i-th argument is an array and the j-entry had an illegal value, then info = -(i*100+j), if the i-th argument is a scalar and had an illegal value, then info = -i.

If info> 0:

If (mod(info,2)0), then one or more eigenvectors failed to converge. Their indices are stored in ifail.

If (mod(info,2,2)0), then eigenvectors corresponding to one or more clusters of eigenvalues could not be reorthogonalized because of insufficient workspace. The indices of the clusters are stored in the array iclustr.

If (mod(info/4,2)0), then space limit prevented p?sygvx from computing all of the eigenvectors between vl and vu. The number of eigenvectors computed is returned in nz.

If (mod(info/8,2)0), then p?stebz failed to compute eigenvalues.

If (mod(info/16,2)0), then B was not positive definite. ifail(1) indicates the order of the smallest minor which is not positive definite.

See Also