Developer Reference for Intel® oneAPI Math Kernel Library for C

ID 766684
Date 11/07/2023
Public

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

Document Table of Contents

?larrv2

Computes the eigenvectors of the tridiagonal matrix T = L*D*LT given L, D and the eigenvalues of L*D*LT.

Syntax

void slarrv2(MKL_INT* n, float* vl, float* vu, float* d, float* l, float* pivmin, MKL_INT* isplit, MKL_INT* m, MKL_INT* dol, MKL_INT* dou, MKL_INT* needil, MKL_INT* neediu, float* minrgp, float* rtol1, float* rtol2, float* w, float* werr, float* wgap, MKL_INT* iblock, MKL_INT* indexw, float* gers, float* sdiam, float* z, MKL_INT* ldz, MKL_INT* isuppz, float* work, MKL_INT* iwork, MKL_INT* vstart, MKL_INT* finish, MKL_INT* maxcls, MKL_INT* ndepth, MKL_INT* parity, MKL_INT* zoffset, MKL_INT* info);

void dlarrv2(MKL_INT* n, double* vl, double* vu, double* d, double* l, double* pivmin, MKL_INT* isplit, MKL_INT* m, MKL_INT* dol, MKL_INT* dou, MKL_INT* needil, MKL_INT* neediu, double* minrgp, double* rtol1, double* rtol2, double* w, double* werr, double* wgap, MKL_INT* iblock, MKL_INT* indexw, double* gers, double* sdiam, double* z, MKL_INT* ldz, MKL_INT* isuppz, double* work, MKL_INT* iwork, MKL_INT* vstart, MKL_INT* finish, MKL_INT* maxcls, MKL_INT* ndepth, MKL_INT* parity, MKL_INT* zoffset, MKL_INT* info);

Include Files

  • mkl_scalapack.h

Description

?larrv2 computes the eigenvectors of the tridiagonal matrix T = LDLT given L, D and approximations to the eigenvalues of LDLT. The input eigenvalues should have been computed by larre2a or by previous calls to ?larrv2.

The major difference between the parallel and the sequential construction of the representation tree is that in the parallel case, not all eigenvalues of a given cluster might be computed locally. Other processors might "own" and refine part of an eigenvalue cluster. This is crucial for scalability. Thus there might be communication necessary before the current level of the representation tree can be parsed.

Please note:

  • The calling sequence has two additional integer parameters, dol and dou, that should satisfy mdoudol1. These parameters are only relevant when both eigenvalues and eigenvectors are computed (stegr2b parameter jobz = 'V'). ?larrv2 only computes the eigenvectors corresponding to eigenvalues dol through dou in w. (That is, instead of computing the eigenvectors belonging to w[0] through w[m-1], only the eigenvectors belonging to eigenvalues w[dol - 1] through w[dou -1] are computed. In this case, only the eigenvalues dol:dou are guaranteed to be accurately refined to all figures by Rayleigh-Quotient iteration.

  • The additional arguments vstart, finish, ndepth, parity, zoffset are included as a thread-safe implementation equivalent to save variables. These variables store details about the local representation tree which is computed layerwise. For scalability reasons, eigenvalues belonging to the locally relevant representation tree might be computed on other processors. These need to be communicated before the inspection of the RRRs can proceed on any given layer. Note that only when the variable finish is non-zero, the computation has ended. All eigenpairs between dol and dou have been computed. m is set to dou - dol + 1.

  • ?larrv2 needs more workspace in z than the sequential slarrv. It is used to store the conformal embedding of the local representation tree.

Product and Performance Information

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

Notice revision #20201201

Input Parameters

n

The order of the matrix. n 0.

vl, vu

Lower and upper bounds of the interval that contains the desired eigenvalues. vl < vu. Needed to compute gaps on the left or right end of the extremal eigenvalues in the desired range. vu is currently not used but kept as parameter in case needed.

d

Array of size n

The n diagonal elements of the diagonal matrix d. On exit, d is overwritten.

l

Array of size n

The (n-1) subdiagonal elements of the unit bidiagonal matrix L are in elements 0 to n-2 of l (if the matrix is not split.) At the end of each block is stored the corresponding shift as given by ?larre. On exit, l is overwritten.

pivmin

The minimum pivot allowed in the sturm sequence.

isplit

Array of size n

The splitting points, at which the matrix T breaks up into blocks. The first block consists of rows/columns 1 to isplit[ 0 ], the second of rows/columns isplit[ 0 ] + 1 through isplit[ 1 ], etc.

m

The total number of input eigenvalues. 0 mn.

dol, dou

If you want to compute only selected eigenvectors from all the eigenvalues supplied, you can specify an index range dol:dou. Or else the setting dol=1, dou=m should be applied. Note that dol and dou refer to the order in which the eigenvalues are stored in w. If you want to compute only selected eigenpairs, the columns dol-1 to dou+1 of the eigenvector space Z contain the computed eigenvectors. All other columns of Z are set to zero.

If dol > 1, then Z(:,dol-1-zoffset) is used.

If dou < m, then Z(:,dou+1-zoffset) is used.

needil, neediu

Describe which are the left and right outermost eigenvalues that still need to be included in the computation. These indices indicate whether eigenvalues from other processors are needed to correctly compute the conformally embedded representation tree.

When dolneedilneediudou, all required eigenvalues are local to the processor and no communication is required to compute its part of the representation tree.

minrgp

The minimum relative gap threshold to decide whether an eigenvalue or a cluster boundary is reached.

rtol1, rtol2

Parameters for bisection. An interval [left,right] has converged if right-left < max( rtol1*gap, rtol2*max(|left|,|right|) )

w

Array of size n

The first m elements of w contain the approximate eigenvalues for which eigenvectors are to be computed. The eigenvalues should be grouped by split-off block and ordered from smallest to largest within the block. (The output array w from ?stegr2a is expected here.) Furthermore, they are with respect to the shift of the corresponding root representation for their block.

werr

Array of size n

The first m elements contain the semiwidth of the uncertainty interval of the corresponding eigenvalue in w.

wgap

Array of size n

The separation from the right neighbor eigenvalue in w.

iblock

Array of size n

The indices of the blocks (submatrices) associated with the corresponding eigenvalues in w; iblock[i]=1 if eigenvalue w[i] belongs to the first block from the top, iblock[i]=2 if w[i] belongs to the second block, and so on.

indexw

Array of size n

The indices of the eigenvalues within each block (submatrix). For example: indexw[i]= 10 and iblock[i]=2 imply that the (i+1)-th eigenvalue w[i] is the 10th eigenvalue in block 2.

gers

Array of size 2*n

The n Gerschgorin intervals (the i-th Gerschgorin interval is (gers[2*i-2], gers[2*i-1])). The Gerschgorin intervals should be computed from the original unshifted matrix.

Not used but kept as parameter for possible future use.

sdiam

Array of size n

The spectral diameters for all unreduced blocks.

ldz

The leading dimension of the array z. ldz 1, and if stegr2b parameter jobz = 'V', ldz max(1,n).

work

(workspace) array of size 12*n

iwork

(workspace) Array of size 7*n

vstart

Non-zero on initialization, set to zero afterwards.

finish

A flag that indicates whether all eigenpairs have been computed.

maxcls

The largest cluster worked on by this processor in the representation tree.

ndepth

The current depth of the representation tree. Set to zero on initial pass, changed when the deeper levels of the representation tree are generated.

parity

An internal parameter needed for the storage of the clusters on the current level of the representation tree.

zoffset

Offset for storing the eigenpairs when z is distributed in 1D-cyclic fashion.

OUTPUT Parameters

needil, neediu

w

Unshifted eigenvalues for which eigenvectors have already been computed.

werr

Contains refined values of its input approximations.

wgap

Contains refined values of its input approximations. Very small gaps are changed.

z

Array of size ldz * max(1,m)

If info = 0, the first m columns of the matrix Z, stored in the array z, contain the orthonormal eigenvectors of the matrix T corresponding to the input eigenvalues, with the i-th column of Z holding the eigenvector associated with w[i - 1].

In the distributed version, only a subset of columns is accessed, see dol, dou, and zoffset.

isuppz

Array of size 2*max(1,m)

The support of the eigenvectors in z, i.e., the indices indicating the non-zero elements in z. The i-th eigenvector is non-zero only in elements isuppz[ 2*i-2 ] through isuppz[ 2*i-1 ].

vstart

Non-zero on initialization, set to zero afterwards.

finish

A flag that indicates whether all eigenpairs have been computed.

maxcls

The largest cluster worked on by this processor in the representation tree.

ndepth

The current depth of the representation tree. Set to zero on initial pass, changed when the deeper levels of the representation tree are generated.

parity

An internal parameter needed for the storage of the clusters on the current level of the representation tree.

info

= 0: successful exit

> 0: A problem occured in ?larrv2.

< 0: One of the called functions signaled an internal problem.

Needs inspection of the corresponding parameter info for further information.

=-1: Problem in ?larrb2 when refining a child's eigenvalues.

=-2: Problem in ?larrf2 when computing the RRR of a child. When a child is inside a tight cluster, it can be difficult to find an RRR. A partial remedy from the user's point of view is to make the parameter minrgp smaller and recompile. However, as the orthogonality of the computed vectors is proportional to 1/minrgp, be aware that decreasing minrgp might be reduce precision.

=-3: Problem in ?larrb2 when refining a single eigenvalue after the Rayleigh correction was rejected.

= 5: The Rayleigh Quotient Iteration failed to converge to full accuracy.

See Also