There is one prototype of syevr
available, please see below.
syevr( const char jobz, const char range, MatrixA& a, const Scalar >, const Scalar >, const int_t il, const int_t iu, const Scalar >, int_t& m, VectorW& w, MatrixZ& z, VectorISUPPZ& isuppz );
syevr (short for $FRIENDLY_NAME)
provides a C++ interface to LAPACK routines SSYEVR and DSYEVR. syevr computes selected eigenvalues
and, optionally, eigenvectors of a real symmetric matrix A. Eigenvalues
and eigenvectors can be selected by specifying either a range of values
or a range of indices for the desired eigenvalues.
syevr first reduces the
matrix A to tridiagonal form T with a call to DSYTRD. Then, whenever
possible, syevr calls
DSTEMR to compute the eigenspectrum using Relatively Robust Representations.
DSTEMR computes eigenvalues by the dqds algorithm, while orthogonal eigenvectors
are computed from various "good" L D L^T representations (also
known as Relatively Robust Representations). Gram-Schmidt orthogonalization
is avoided as far as possible. More specifically, the various steps of
the algorithm are as follows.
For each unreduced block (submatrix) of T, (a) Compute T - sigma I = L D L^T, so that L and D define all the wanted eigenvalues to high relative accuracy. This means that small relative changes in the entries of D and L cause only small relative changes in the eigenvalues and eigenvectors. The standard (unfactored) representation of the tridiagonal matrix T does not have this property in general. (b) Compute the eigenvalues to suitable accuracy. If the eigenvectors are desired, the algorithm attains full accuracy of the computed eigenvalues only right before the corresponding vectors have to be computed, see steps c) and d). (c) For each cluster of close eigenvalues, select a new shift close to the cluster, find a new factorization, and refine the shifted eigenvalues to suitable accuracy. (d) For each eigenvalue with a large enough relative separation compute the corresponding eigenvector by forming a rank revealing twisted factorization. Go back to (c) for any clusters that remain.
The desired accuracy of the output can be specified by the input parameter ABSTOL.
For more details, see DSTEMR's documentation and: - Inderjit S. Dhillon and Beresford N. Parlett: "Multiple representations to compute orthogonal eigenvectors of symmetric tridiagonal matrices," Linear Algebra and its Applications, 387(1), pp. 1-28, August 2004. - Inderjit Dhillon and Beresford Parlett: "Orthogonal Eigenvectors and Relative Gaps," SIAM Journal on Matrix Analysis and Applications, Vol. 25, 2004. Also LAPACK Working Note 154. - Inderjit Dhillon: "A new O(n^2) algorithm for the symmetric tridiagonal eigenvalue/eigenvector problem", Computer Science Division Technical Report No. UCB/CSD-97-971, UC Berkeley, May 1997.
Note 1 : syevr calls
DSTEMR when the full spectrum is requested on machines which conform
to the ieee-754 floating point standard. syevr
calls DSTEBZ and SSTEIN on non-ieee machines and when partial spectrum
requests are made.
Normal execution of DSTEMR may create NaNs and infinities and hence may abort due to a floating point exception in environments which do not handle NaNs and infinities in the ieee standard default manner.
The selection of the LAPACK routine is done during compile-time, and
is determined by the type of values contained in type MatrixA.
The type of values is obtained through the value_type
meta-function typename value_type<MatrixA>::type. The dispatching table below illustrates
to which specific routine the code path will be generated.
Defined in header boost/numeric/bindings/lapack/driver/syevr.hpp.
Parameters
The definition of term 1
The definition of term 2
The definition of term 3.
Definitions may contain paragraphs.
#include <boost/numeric/bindings/lapack/driver/syevr.hpp> using namespace boost::numeric::bindings; lapack::syevr( x, y, z );
this will output
[5] 0 1 2 3 4 5