Module: Linear Algebra Group: Symmetric Eigenvalue Decomposition classes
Does not inherit
#include <rw/lapack/seigsrv.h> RWSymSomeEigServer<double> server; RWSymEigDecomp<double> eig = server(A); // A is a
// RWSymBandMat<double> #include <rw/lapack/heigsrv.h> RWHermSomeEigServer<double> server; RWHermEigDecomp<double> eig = server(A); // A is a
// RWHermBandMat<DComplex>
The symmetric eigenvalue server class, RWSymSomeEigServer<T>, and the Hermitian eigenvalue server class, RWHermSomeEigServer<T>, allow the computation of a subset of the eigenvalues and (optionally) their corresponding eigenvectors. The computation uses the bisection method to find the eigenvalues, followed by inverse iteration to find the eigenvectors. The subset of eigenvalues to be computed is specified using the RWSlice class, or one of its subclasses, RWRange or RWToEnd. This provides the flexibility to specify either the smallest eigenvalues, the largest, or a selection in between. The eigenvalue ordering is smallest to largest.
#include <iostream> #include <rw/lapack/seigsrv.h> int main() { RWSymMat<double> A; // input a matrix cin >> A; DoubleSomeEigServer server; // configure a server server.setRange(RWSlice(0,5)); // the 5 smallest
eigenvalues RWSymEigDecomp<double> eig = server(A); return 0; }
RWSymSomeEigServer(bool computeVecs=true); RWHermSomeEigServer(bool computeVecs=true);
Constructs a server. The parameter indicates whether this server should be configured to compute eigenvectors as well as eigenvalues.
void RWSymSomeEigServer<T>::computeEigenVectors (bool); void RWHermSomeEigServer<DComplex>::computeEigenVectors
(bool);
Sets whether or not the server should compute eigenvectors as well as eigenvalues.
bool RWSymSomeEigServer<T>::computeEigenVectors() const; bool RWHermSomeEigServer<DComplex>::computeEigenVectors() const;
Returns true if this server is configured to compute eigenvectors as well as eigenvalues.
RWSymEigDecomp<T> RWSymSomeEigServer<T>::decompose
(const FloatSymTriDiagDecomp& A) const; DoubleSymEigDecomp<DComplex> RWHermSomeEigServer<DComplex>::decompose
(const DoubleSymTriDiagDecomp& A) const;
Computes the eigenvalue decomposition of the tridiagonal matrix inside the tridiagonal decomposition.
RWSlice RWSymSomeEigServer<T>::setRange(const RWSlice&); RWSlice RWHermSomeEigServer<T>::setRange(const RWSlice&);
Sets the range of eigenvalues to be computed. Returns the previous range.
T RWSymSomeEigServer<T>::setTolerance(T); double RWHermSomeEigServer<DComplex>::setTolerance(double);
Sets the tolerance to which we have to compute the eigenvalues. Returns the previous tolerance.
RWSymEigDecomp<T> RWSymSomeEigServer<T>::operator()
(const RWSymMat<T>& A) const; RWSymEigDecomp<T> RWSymSomeEigServer<T>::operator()
(const RWSymBandMat<T>& A) const; RWHermEigDecomp<DComplex> RWHermSomeEigServer<DComplex>::operator()
(const RWHermMat<DComplex>& A) const; RWHermEigDecomp<DComplex> RWHermSomeEigServer<DComplex>::operator()
(const RWHermBandMat<Dcomplex>& A) const;
Computes a symmetric eigenvalue decomposition.
© Copyright Rogue Wave Software, Inc. All Rights Reserved.
Rogue Wave and SourcePro are registered trademarks of Rogue Wave Software, Inc. in the United States and other countries. All other trademarks are the property of their respective owners.
Contact Rogue Wave about documentation or support issues.