>
>#include <rw/fseigsrv.h> // FloatSymSomeEigServer #include <rw/dseigsrv.h> // DoubleSymSomeEigServer #include <rw/cheigsrv.h> // DComplexHermSomeEigServer DoubleSymSomeEigServer server; DoubleSymEigDecomp eig = server(A); // A is a
// DoubleSym[Band]Mat
The symmetric eigenvalue server classes, {TYPE}SymSomeEigServer, and the Hermitian eigenvalue server class, DComplexHermSomeEigServer, 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 <rw/dseigsrv.h> main() { DoubleSymMat A; // input a matrix cin >> A; DoubleSomeEigServer server; // configure a server server.setRange(RWSlice(0,5)); // the 5 smallest eigenvalues DoubleSymEigDecomp eig = server(A); }>
FloatSymSomeEigServer(RWBoolean computeVecs=TRUE); DoubleSymSomeEigServer(RWBoolean computeVecs=TRUE); DComplexHermSomeEigServer(RWBoolean computeVecs=TRUE);
Constructs a server. The parameter indicates whether this server should be configured to compute eigenvectors as well as eigenvalues.
void FloatSymSomeEigServer::computeEigenVectors(RWBoolean); void DoubleSymSomeEigServer::computeEigenVectors(RWBoolean); void DComplexHermSomeEigServer::computeEigenVectors(RWBoolean);
Sets whether or not the server should compute eigenvectors as well as eigenvalues.
RWBoolean FloatSymSomeEigServer::computeEigenVectors() const; RWBoolean DoubleSymSomeEigServer::computeEigenVectors() const; RWBoolean DComplexHermSomeEigServer::computeEigenVectors() const;
Returns TRUE if this server is configured to compute eigenvectors as well as eigenvalues.
FloatSymEigDecomp FloatSymSomeEigServer::decompose
(const FloatSymTriDiagDecomp& A) const; DoubleSymEigDecomp DoubleSymSomeEigServer::decompose
(const DoubleSymTriDiagDecomp& A) const; DoubleSymEigDecomp DComplexHermSomeEigServer::decompose
(const DoubleSymTriDiagDecomp& A) const;
Computes the eigenvalue decomposition of the tridiagonal matrix inside the tridiagonal decomposition.
RWSlice FloatSymSomeEigServer::setRange(const RWSlice&); RWSlice DoubleSymSomeEigServer::setRange(const RWSlice&); RWSlice DComplexHermSomeEigServer::setRange(const RWSlice&);
Sets the range of eigenvalues to be computed. Returns the previous range.
float FloatSymSomeEigServer::setTolerance(float); double DoubleSymSomeEigServer::setTolerance(double); double DComplexHermSomeEigServer::setTolerance(double);
Sets the tolerance to which we have to compute the eigenvalues. Returns the previous tolerance.
FloatSymEigDecomp FloatSymSomeEigDecomp::operator()
(const FloatSymMat& A) const; FloatSymEigDecomp FloatSymSomeEigDecomp::operator()
(const FloatSymBandMat& A) const; DoubleSymEigDecomp DoubleSymSomeEigDecomp::operator()
(const DoubleSymMat& A) const; DoubleSymEigDecomp DoubleSymSomeEigDecomp::operator()
(const DoubleSymBandMat& A) const; DComplexSymEigDecomp DComplexSymSomeEigDecomp::operator()
(const DComplexSymMat& A) const; DComplexHermEigDecomp DComplexHermSomeEigDecomp::operator()
(const DComplexHermBandMat& A) const;
Computes a symmetric eigenvalue decomposition.
©Copyright 1999, Rogue Wave Software, Inc.
Send mail to report errors or comment on the documentation.