10#ifndef EIGEN_ARPACKGENERALIZEDSELFADJOINTEIGENSOLVER_H
11#define EIGEN_ARPACKGENERALIZEDSELFADJOINTEIGENSOLVER_H
13#include "../../../../Eigen/Dense"
16#include "./InternalHeaderCheck.h"
21template <
typename Scalar,
typename RealScalar>
23template <
typename MatrixSolver,
typename MatrixType,
typename Scalar,
bool BisSPD>
27template <
typename MatrixType,
typename MatrixSolver = SimplicialLLT<MatrixType>,
bool BisSPD = false>
28class ArpackGeneralizedSelfAdjointEigenSolver {
33 typedef typename MatrixType::Scalar Scalar;
34 typedef typename MatrixType::Index Index;
42 typedef typename NumTraits<Scalar>::Real RealScalar;
49 typedef typename internal::plain_col_type<MatrixType, RealScalar>::type RealVectorType;
57 ArpackGeneralizedSelfAdjointEigenSolver()
60 m_isInitialized(false),
61 m_eigenvectorsOk(false),
87 ArpackGeneralizedSelfAdjointEigenSolver(
const MatrixType &A,
const MatrixType &B, Index nbrEigenvalues,
92 m_isInitialized(false),
93 m_eigenvectorsOk(false),
96 compute(A, B, nbrEigenvalues, eigs_sigma, options, tol);
121 ArpackGeneralizedSelfAdjointEigenSolver(
const MatrixType &A, Index nbrEigenvalues, std::string eigs_sigma =
"LM",
125 m_isInitialized(false),
126 m_eigenvectorsOk(false),
129 compute(A, nbrEigenvalues, eigs_sigma, options, tol);
155 ArpackGeneralizedSelfAdjointEigenSolver &compute(
const MatrixType &A,
const MatrixType &B, Index nbrEigenvalues,
157 RealScalar tol = 0.0);
181 ArpackGeneralizedSelfAdjointEigenSolver &compute(
const MatrixType &A, Index nbrEigenvalues,
183 RealScalar tol = 0.0);
204 const Matrix<Scalar, Dynamic, Dynamic> &eigenvectors()
const {
205 eigen_assert(m_isInitialized &&
"ArpackGeneralizedSelfAdjointEigenSolver is not initialized.");
206 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
225 const Matrix<Scalar, Dynamic, 1> &eigenvalues()
const {
226 eigen_assert(m_isInitialized &&
"ArpackGeneralizedSelfAdjointEigenSolver is not initialized.");
248 Matrix<Scalar, Dynamic, Dynamic> operatorSqrt()
const {
249 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
250 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
251 return m_eivec * m_eivalues.cwiseSqrt().asDiagonal() * m_eivec.adjoint();
272 Matrix<Scalar, Dynamic, Dynamic> operatorInverseSqrt()
const {
273 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
274 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
275 return m_eivec * m_eivalues.cwiseInverse().cwiseSqrt().asDiagonal() * m_eivec.adjoint();
283 eigen_assert(m_isInitialized &&
"ArpackGeneralizedSelfAdjointEigenSolver is not initialized.");
287 size_t getNbrConvergedEigenValues()
const {
return m_nbrConverged; }
289 size_t getNbrIterations()
const {
return m_nbrIterations; }
292 Matrix<Scalar, Dynamic, Dynamic> m_eivec;
293 Matrix<Scalar, Dynamic, 1> m_eivalues;
295 bool m_isInitialized;
296 bool m_eigenvectorsOk;
298 size_t m_nbrConverged;
299 size_t m_nbrIterations;
302template <
typename MatrixType,
typename MatrixSolver,
bool BisSPD>
303ArpackGeneralizedSelfAdjointEigenSolver<MatrixType, MatrixSolver, BisSPD> &
304ArpackGeneralizedSelfAdjointEigenSolver<MatrixType, MatrixSolver, BisSPD>::compute(
const MatrixType &A,
305 Index nbrEigenvalues,
306 std::string eigs_sigma,
int options,
309 compute(A, B, nbrEigenvalues, eigs_sigma, options, tol);
314template <
typename MatrixType,
typename MatrixSolver,
bool BisSPD>
315ArpackGeneralizedSelfAdjointEigenSolver<MatrixType, MatrixSolver, BisSPD> &
316ArpackGeneralizedSelfAdjointEigenSolver<MatrixType, MatrixSolver, BisSPD>::compute(
const MatrixType &A,
318 Index nbrEigenvalues,
319 std::string eigs_sigma,
int options,
321 eigen_assert(A.cols() == A.rows());
322 eigen_assert(B.cols() == B.rows());
323 eigen_assert(B.rows() == 0 || A.cols() == B.rows());
324 eigen_assert((options & ~(EigVecMask | GenEigMask)) == 0 && (options & EigVecMask) != EigVecMask &&
325 "invalid option parameter");
327 bool isBempty = (B.rows() == 0) || (B.cols() == 0);
335 int n = (int)A.cols();
343 RealScalar sigma = 0.0;
345 if (eigs_sigma.length() >= 2 && isalpha(eigs_sigma[0]) && isalpha(eigs_sigma[1])) {
346 eigs_sigma[0] = toupper(eigs_sigma[0]);
347 eigs_sigma[1] = toupper(eigs_sigma[1]);
353 if (eigs_sigma.substr(0, 2) !=
"SM") {
354 whch[0] = eigs_sigma[0];
355 whch[1] = eigs_sigma[1];
358 eigen_assert(
false &&
"Specifying clustered eigenvalues is not yet supported!");
363 sigma = atof(eigs_sigma.c_str());
372 if (eigs_sigma.substr(0, 2) ==
"SM" || !(isalpha(eigs_sigma[0]) && isalpha(eigs_sigma[1])) || (!isBempty && !BisSPD))
377 int mode = (bmat[0] ==
'G') + 1;
378 if (eigs_sigma.substr(0, 2) ==
"SM" || !(isalpha(eigs_sigma[0]) && isalpha(eigs_sigma[1]))) {
387 int nev = (int)nbrEigenvalues;
397 int ncv = std::min(std::max(2 * nev, 20), n);
407 int lworkl = ncv * ncv + 8 * ncv;
410 int *iparam =
new int[11];
412 iparam[2] = std::max(300, (
int)std::ceil(2 * n / std::max(ncv, 1)));
417 int *ipntr =
new int[11];
437 if (mode == 1 || mode == 2) {
438 if (!isBempty) OP.compute(B);
439 }
else if (mode == 3) {
446 MatrixType AminusSigmaB(A);
447 for (
Index i = 0; i < A.rows(); ++i) AminusSigmaB.coeffRef(i, i) -= sigma;
449 OP.compute(AminusSigmaB);
451 MatrixType AminusSigmaB = A - sigma * B;
452 OP.compute(AminusSigmaB);
457 if (!(mode == 1 && isBempty) && !(mode == 2 && isBempty) && OP.info() !=
Success)
458 std::cout <<
"Error factoring matrix" << std::endl;
461 internal::arpack_wrapper<Scalar, RealScalar>::saupd(&ido, bmat, &n, whch, &nev, &tol, resid, &ncv, v, &ldv, iparam,
462 ipntr, workd, workl, &lworkl, &info);
464 if (ido == -1 || ido == 1) {
465 Scalar *in = workd + ipntr[0] - 1;
466 Scalar *out = workd + ipntr[1] - 1;
468 if (ido == 1 && mode != 2) {
469 Scalar *out2 = workd + ipntr[2] - 1;
470 if (isBempty || mode == 1)
475 in = workd + ipntr[2] - 1;
486 internal::OP<MatrixSolver, MatrixType, Scalar, BisSPD>::applyOP(OP, A, n, in, out);
488 }
else if (mode == 2) {
494 }
else if (mode == 3) {
498 if (ido == 1 || isBempty)
503 }
else if (ido == 2) {
504 Scalar *in = workd + ipntr[0] - 1;
505 Scalar *out = workd + ipntr[1] - 1;
507 if (isBempty || mode == 1)
521 eigen_assert(
false &&
"Unknown ARPACK return value!");
529 char howmny[2] =
"A";
533 int *select =
new int[ncv];
537 m_eivalues.resize(nev, 1);
539 internal::arpack_wrapper<Scalar, RealScalar>::seupd(&rvec, howmny, select, m_eivalues.data(), v, &ldv, &sigma, bmat,
540 &n, whch, &nev, &tol, resid, &ncv, v, &ldv, iparam, ipntr,
541 workd, workl, &lworkl, &info);
549 m_eivec.resize(A.rows(), nev);
550 for (
int i = 0; i < nev; i++)
551 for (
int j = 0; j < n; j++) m_eivec(j, i) = v[i * n + j] / scale;
553 if (mode == 1 && !isBempty && BisSPD)
554 internal::OP<MatrixSolver, MatrixType, Scalar, BisSPD>::project(OP, n, nev, m_eivec.data());
556 m_eigenvectorsOk =
true;
559 m_nbrIterations = iparam[2];
560 m_nbrConverged = iparam[4];
575 m_isInitialized =
true;
582extern "C" void ssaupd_(
int *ido,
char *bmat,
int *n,
char *which,
int *nev,
float *tol,
float *resid,
int *ncv,
583 float *v,
int *ldv,
int *iparam,
int *ipntr,
float *workd,
float *workl,
int *lworkl,
586extern "C" void sseupd_(
int *rvec,
char *All,
int *select,
float *d,
float *z,
int *ldz,
float *sigma,
char *bmat,
587 int *n,
char *which,
int *nev,
float *tol,
float *resid,
int *ncv,
float *v,
int *ldv,
588 int *iparam,
int *ipntr,
float *workd,
float *workl,
int *lworkl,
int *ierr);
592extern "C" void dsaupd_(
int *ido,
char *bmat,
int *n,
char *which,
int *nev,
double *tol,
double *resid,
int *ncv,
593 double *v,
int *ldv,
int *iparam,
int *ipntr,
double *workd,
double *workl,
int *lworkl,
596extern "C" void dseupd_(
int *rvec,
char *All,
int *select,
double *d,
double *z,
int *ldz,
double *sigma,
char *bmat,
597 int *n,
char *which,
int *nev,
double *tol,
double *resid,
int *ncv,
double *v,
int *ldv,
598 int *iparam,
int *ipntr,
double *workd,
double *workl,
int *lworkl,
int *ierr);
602template <
typename Scalar,
typename RealScalar>
603struct arpack_wrapper {
604 static inline void saupd(
int *ido,
char *bmat,
int *n,
char *which,
int *nev, RealScalar *tol, Scalar *resid,
605 int *ncv, Scalar *v,
int *ldv,
int *iparam,
int *ipntr, Scalar *workd, Scalar *workl,
606 int *lworkl,
int *info) {
607 EIGEN_STATIC_ASSERT(!NumTraits<Scalar>::IsComplex, NUMERIC_TYPE_MUST_BE_REAL)
610 static inline void seupd(
int *rvec,
char *All,
int *select, Scalar *d, Scalar *z,
int *ldz, RealScalar *sigma,
611 char *bmat,
int *n,
char *which,
int *nev, RealScalar *tol, Scalar *resid,
int *ncv,
612 Scalar *v,
int *ldv,
int *iparam,
int *ipntr, Scalar *workd, Scalar *workl,
int *lworkl,
614 EIGEN_STATIC_ASSERT(!NumTraits<Scalar>::IsComplex, NUMERIC_TYPE_MUST_BE_REAL)
619struct arpack_wrapper<float, float> {
620 static inline void saupd(
int *ido,
char *bmat,
int *n,
char *which,
int *nev,
float *tol,
float *resid,
int *ncv,
621 float *v,
int *ldv,
int *iparam,
int *ipntr,
float *workd,
float *workl,
int *lworkl,
623 ssaupd_(ido, bmat, n, which, nev, tol, resid, ncv, v, ldv, iparam, ipntr, workd, workl, lworkl, info);
626 static inline void seupd(
int *rvec,
char *All,
int *select,
float *d,
float *z,
int *ldz,
float *sigma,
char *bmat,
627 int *n,
char *which,
int *nev,
float *tol,
float *resid,
int *ncv,
float *v,
int *ldv,
628 int *iparam,
int *ipntr,
float *workd,
float *workl,
int *lworkl,
int *ierr) {
629 sseupd_(rvec, All, select, d, z, ldz, sigma, bmat, n, which, nev, tol, resid, ncv, v, ldv, iparam, ipntr, workd,
630 workl, lworkl, ierr);
635struct arpack_wrapper<double, double> {
636 static inline void saupd(
int *ido,
char *bmat,
int *n,
char *which,
int *nev,
double *tol,
double *resid,
int *ncv,
637 double *v,
int *ldv,
int *iparam,
int *ipntr,
double *workd,
double *workl,
int *lworkl,
639 dsaupd_(ido, bmat, n, which, nev, tol, resid, ncv, v, ldv, iparam, ipntr, workd, workl, lworkl, info);
642 static inline void seupd(
int *rvec,
char *All,
int *select,
double *d,
double *z,
int *ldz,
double *sigma,
char *bmat,
643 int *n,
char *which,
int *nev,
double *tol,
double *resid,
int *ncv,
double *v,
int *ldv,
644 int *iparam,
int *ipntr,
double *workd,
double *workl,
int *lworkl,
int *ierr) {
645 dseupd_(rvec, All, select, d, v, ldv, sigma, bmat, n, which, nev, tol, resid, ncv, v, ldv, iparam, ipntr, workd,
646 workl, lworkl, ierr);
650template <
typename MatrixSolver,
typename MatrixType,
typename Scalar,
bool BisSPD>
652 static inline void applyOP(MatrixSolver &OP,
const MatrixType &A,
int n, Scalar *in, Scalar *out);
653 static inline void project(MatrixSolver &OP,
int n,
int k, Scalar *vecs);
656template <
typename MatrixSolver,
typename MatrixType,
typename Scalar>
657struct OP<MatrixSolver, MatrixType, Scalar, true> {
658 static inline void applyOP(MatrixSolver &OP,
const MatrixType &A,
int n, Scalar *in, Scalar *out) {
663 Matrix<Scalar, Dynamic, 1>::Map(out, n) = OP.matrixU().solve(Matrix<Scalar, Dynamic, 1>::Map(in, n));
664 Matrix<Scalar, Dynamic, 1>::Map(out, n) = OP.permutationPinv() * Matrix<Scalar, Dynamic, 1>::Map(out, n);
668 Matrix<Scalar, Dynamic, 1>::Map(out, n) = A * Matrix<Scalar, Dynamic, 1>::Map(out, n);
672 Matrix<Scalar, Dynamic, 1>::Map(out, n) = OP.permutationP() * Matrix<Scalar, Dynamic, 1>::Map(out, n);
673 Matrix<Scalar, Dynamic, 1>::Map(out, n) = OP.matrixL().solve(Matrix<Scalar, Dynamic, 1>::Map(out, n));
676 static inline void project(MatrixSolver &OP,
int n,
int k, Scalar *vecs) {
679 Matrix<Scalar, Dynamic, Dynamic>::Map(vecs, n, k) =
680 OP.matrixU().solve(Matrix<Scalar, Dynamic, Dynamic>::Map(vecs, n, k));
681 Matrix<Scalar, Dynamic, Dynamic>::Map(vecs, n, k) =
682 OP.permutationPinv() * Matrix<Scalar, Dynamic, Dynamic>::Map(vecs, n, k);
686template <
typename MatrixSolver,
typename MatrixType,
typename Scalar>
687struct OP<MatrixSolver, MatrixType, Scalar, false> {
688 static inline void applyOP(MatrixSolver &OP,
const MatrixType &A,
int n, Scalar *in, Scalar *out) {
689 eigen_assert(
false &&
"Should never be in here...");
692 static inline void project(MatrixSolver &OP,
int n,
int k, Scalar *vecs) {
693 eigen_assert(
false &&
"Should never be in here...");
Namespace containing all symbols from the Eigen library.
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index