Eigen-unsupported  5.0.0-dev
Loading...
Searching...
No Matches
ArpackSelfAdjointEigenSolver.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2012 David Harmon <dharmon@gmail.com>
5//
6// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10#ifndef EIGEN_ARPACKGENERALIZEDSELFADJOINTEIGENSOLVER_H
11#define EIGEN_ARPACKGENERALIZEDSELFADJOINTEIGENSOLVER_H
12
13#include "../../../../Eigen/Dense"
14
15// IWYU pragma: private
16#include "./InternalHeaderCheck.h"
17
18namespace Eigen {
19
20namespace internal {
21template <typename Scalar, typename RealScalar>
22struct arpack_wrapper;
23template <typename MatrixSolver, typename MatrixType, typename Scalar, bool BisSPD>
24struct OP;
25} // namespace internal
26
27template <typename MatrixType, typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
28class ArpackGeneralizedSelfAdjointEigenSolver {
29 public:
30 // typedef typename MatrixSolver::MatrixType MatrixType;
31
33 typedef typename MatrixType::Scalar Scalar;
34 typedef typename MatrixType::Index Index;
35
42 typedef typename NumTraits<Scalar>::Real RealScalar;
43
49 typedef typename internal::plain_col_type<MatrixType, RealScalar>::type RealVectorType;
50
57 ArpackGeneralizedSelfAdjointEigenSolver()
58 : m_eivec(),
59 m_eivalues(),
60 m_isInitialized(false),
61 m_eigenvectorsOk(false),
62 m_nbrConverged(0),
63 m_nbrIterations(0) {}
64
87 ArpackGeneralizedSelfAdjointEigenSolver(const MatrixType &A, const MatrixType &B, Index nbrEigenvalues,
88 std::string eigs_sigma = "LM", int options = ComputeEigenvectors,
89 RealScalar tol = 0.0)
90 : m_eivec(),
91 m_eivalues(),
92 m_isInitialized(false),
93 m_eigenvectorsOk(false),
94 m_nbrConverged(0),
95 m_nbrIterations(0) {
96 compute(A, B, nbrEigenvalues, eigs_sigma, options, tol);
97 }
98
120
121 ArpackGeneralizedSelfAdjointEigenSolver(const MatrixType &A, Index nbrEigenvalues, std::string eigs_sigma = "LM",
122 int options = ComputeEigenvectors, RealScalar tol = 0.0)
123 : m_eivec(),
124 m_eivalues(),
125 m_isInitialized(false),
126 m_eigenvectorsOk(false),
127 m_nbrConverged(0),
128 m_nbrIterations(0) {
129 compute(A, nbrEigenvalues, eigs_sigma, options, tol);
130 }
131
155 ArpackGeneralizedSelfAdjointEigenSolver &compute(const MatrixType &A, const MatrixType &B, Index nbrEigenvalues,
156 std::string eigs_sigma = "LM", int options = ComputeEigenvectors,
157 RealScalar tol = 0.0);
158
181 ArpackGeneralizedSelfAdjointEigenSolver &compute(const MatrixType &A, Index nbrEigenvalues,
182 std::string eigs_sigma = "LM", int options = ComputeEigenvectors,
183 RealScalar tol = 0.0);
184
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.");
207 return m_eivec;
208 }
209
225 const Matrix<Scalar, Dynamic, 1> &eigenvalues() const {
226 eigen_assert(m_isInitialized && "ArpackGeneralizedSelfAdjointEigenSolver is not initialized.");
227 return m_eivalues;
228 }
229
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();
252 }
253
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();
276 }
277
282 ComputationInfo info() const {
283 eigen_assert(m_isInitialized && "ArpackGeneralizedSelfAdjointEigenSolver is not initialized.");
284 return m_info;
285 }
286
287 size_t getNbrConvergedEigenValues() const { return m_nbrConverged; }
288
289 size_t getNbrIterations() const { return m_nbrIterations; }
290
291 protected:
292 Matrix<Scalar, Dynamic, Dynamic> m_eivec;
293 Matrix<Scalar, Dynamic, 1> m_eivalues;
294 ComputationInfo m_info;
295 bool m_isInitialized;
296 bool m_eigenvectorsOk;
297
298 size_t m_nbrConverged;
299 size_t m_nbrIterations;
300};
301
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,
307 RealScalar tol) {
308 MatrixType B(0, 0);
309 compute(A, B, nbrEigenvalues, eigs_sigma, options, tol);
310
311 return *this;
312}
313
314template <typename MatrixType, typename MatrixSolver, bool BisSPD>
315ArpackGeneralizedSelfAdjointEigenSolver<MatrixType, MatrixSolver, BisSPD> &
316ArpackGeneralizedSelfAdjointEigenSolver<MatrixType, MatrixSolver, BisSPD>::compute(const MatrixType &A,
317 const MatrixType &B,
318 Index nbrEigenvalues,
319 std::string eigs_sigma, int options,
320 RealScalar tol) {
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");
326
327 bool isBempty = (B.rows() == 0) || (B.cols() == 0);
328
329 // For clarity, all parameters match their ARPACK name
330 //
331 // Always 0 on the first call
332 //
333 int ido = 0;
334
335 int n = (int)A.cols();
336
337 // User options: "LA", "SA", "SM", "LM", "BE"
338 //
339 char whch[3] = "LM";
340
341 // Specifies the shift if iparam[6] = { 3, 4, 5 }, not used if iparam[6] = { 1, 2 }
342 //
343 RealScalar sigma = 0.0;
344
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]);
348
349 // In the following special case we're going to invert the problem, since solving
350 // for larger magnitude is much much faster
351 // i.e., if 'SM' is specified, we're going to really use 'LM', the default
352 //
353 if (eigs_sigma.substr(0, 2) != "SM") {
354 whch[0] = eigs_sigma[0];
355 whch[1] = eigs_sigma[1];
356 }
357 } else {
358 eigen_assert(false && "Specifying clustered eigenvalues is not yet supported!");
359
360 // If it's not scalar values, then the user may be explicitly
361 // specifying the sigma value to cluster the evs around
362 //
363 sigma = atof(eigs_sigma.c_str());
364
365 // If atof fails, it returns 0.0, which is a fine default
366 //
367 }
368
369 // "I" means normal eigenvalue problem, "G" means generalized
370 //
371 char bmat[2] = "I";
372 if (eigs_sigma.substr(0, 2) == "SM" || !(isalpha(eigs_sigma[0]) && isalpha(eigs_sigma[1])) || (!isBempty && !BisSPD))
373 bmat[0] = 'G';
374
375 // Now we determine the mode to use
376 //
377 int mode = (bmat[0] == 'G') + 1;
378 if (eigs_sigma.substr(0, 2) == "SM" || !(isalpha(eigs_sigma[0]) && isalpha(eigs_sigma[1]))) {
379 // We're going to use shift-and-invert mode, and basically find
380 // the largest eigenvalues of the inverse operator
381 //
382 mode = 3;
383 }
384
385 // The user-specified number of eigenvalues/vectors to compute
386 //
387 int nev = (int)nbrEigenvalues;
388
389 // Allocate space for ARPACK to store the residual
390 //
391 Scalar *resid = new Scalar[n];
392
393 // Number of Lanczos vectors, must satisfy nev < ncv <= n
394 // Note that this indicates that nev != n, and we cannot compute
395 // all eigenvalues of a mtrix
396 //
397 int ncv = std::min(std::max(2 * nev, 20), n);
398
399 // The working n x ncv matrix, also store the final eigenvectors (if computed)
400 //
401 Scalar *v = new Scalar[n * ncv];
402 int ldv = n;
403
404 // Working space
405 //
406 Scalar *workd = new Scalar[3 * n];
407 int lworkl = ncv * ncv + 8 * ncv; // Must be at least this length
408 Scalar *workl = new Scalar[lworkl];
409
410 int *iparam = new int[11];
411 iparam[0] = 1; // 1 means we let ARPACK perform the shifts, 0 means we'd have to do it
412 iparam[2] = std::max(300, (int)std::ceil(2 * n / std::max(ncv, 1)));
413 iparam[6] = mode; // The mode, 1 is standard ev problem, 2 for generalized ev, 3 for shift-and-invert
414
415 // Used during reverse communicate to notify where arrays start
416 //
417 int *ipntr = new int[11];
418
419 // Error codes are returned in here, initial value of 0 indicates a random initial
420 // residual vector is used, any other values means resid contains the initial residual
421 // vector, possibly from a previous run
422 //
423 int info = 0;
424
425 Scalar scale = 1.0;
426 // if (!isBempty)
427 //{
428 // Scalar scale = B.norm() / std::sqrt(n);
429 // scale = std::pow(2, std::floor(std::log(scale+1)));
431 // for (size_t i=0; i<(size_t)B.outerSize(); i++)
432 // for (typename MatrixType::InnerIterator it(B, i); it; ++it)
433 // it.valueRef() /= scale;
434 // }
435
436 MatrixSolver OP;
437 if (mode == 1 || mode == 2) {
438 if (!isBempty) OP.compute(B);
439 } else if (mode == 3) {
440 if (sigma == 0.0) {
441 OP.compute(A);
442 } else {
443 // Note: We will never enter here because sigma must be 0.0
444 //
445 if (isBempty) {
446 MatrixType AminusSigmaB(A);
447 for (Index i = 0; i < A.rows(); ++i) AminusSigmaB.coeffRef(i, i) -= sigma;
448
449 OP.compute(AminusSigmaB);
450 } else {
451 MatrixType AminusSigmaB = A - sigma * B;
452 OP.compute(AminusSigmaB);
453 }
454 }
455 }
456
457 if (!(mode == 1 && isBempty) && !(mode == 2 && isBempty) && OP.info() != Success)
458 std::cout << "Error factoring matrix" << std::endl;
459
460 do {
461 internal::arpack_wrapper<Scalar, RealScalar>::saupd(&ido, bmat, &n, whch, &nev, &tol, resid, &ncv, v, &ldv, iparam,
462 ipntr, workd, workl, &lworkl, &info);
463
464 if (ido == -1 || ido == 1) {
465 Scalar *in = workd + ipntr[0] - 1;
466 Scalar *out = workd + ipntr[1] - 1;
467
468 if (ido == 1 && mode != 2) {
469 Scalar *out2 = workd + ipntr[2] - 1;
470 if (isBempty || mode == 1)
472 else
474
475 in = workd + ipntr[2] - 1;
476 }
477
478 if (mode == 1) {
479 if (isBempty) {
480 // OP = A
481 //
483 } else {
484 // OP = L^{-1}AL^{-T}
485 //
486 internal::OP<MatrixSolver, MatrixType, Scalar, BisSPD>::applyOP(OP, A, n, in, out);
487 }
488 } else if (mode == 2) {
490
491 // OP = B^{-1} A
492 //
494 } else if (mode == 3) {
495 // OP = (A-\sigmaB)B (\sigma could be 0, and B could be I)
496 // The B * in is already computed and stored at in if ido == 1
497 //
498 if (ido == 1 || isBempty)
500 else
502 }
503 } else if (ido == 2) {
504 Scalar *in = workd + ipntr[0] - 1;
505 Scalar *out = workd + ipntr[1] - 1;
506
507 if (isBempty || mode == 1)
509 else
511 }
512 } while (ido != 99);
513
514 if (info == 1)
515 m_info = NoConvergence;
516 else if (info == 3)
517 m_info = NumericalIssue;
518 else if (info < 0)
519 m_info = InvalidInput;
520 else if (info != 0)
521 eigen_assert(false && "Unknown ARPACK return value!");
522 else {
523 // Do we compute eigenvectors or not?
524 //
525 int rvec = (options & ComputeEigenvectors) == ComputeEigenvectors;
526
527 // "A" means "All", use "S" to choose specific eigenvalues (not yet supported in ARPACK))
528 //
529 char howmny[2] = "A";
530
531 // if howmny == "S", specifies the eigenvalues to compute (not implemented in ARPACK)
532 //
533 int *select = new int[ncv];
534
535 // Final eigenvalues
536 //
537 m_eivalues.resize(nev, 1);
538
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);
542
543 if (info == -14)
544 m_info = NoConvergence;
545 else if (info != 0)
546 m_info = InvalidInput;
547 else {
548 if (rvec) {
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;
552
553 if (mode == 1 && !isBempty && BisSPD)
554 internal::OP<MatrixSolver, MatrixType, Scalar, BisSPD>::project(OP, n, nev, m_eivec.data());
555
556 m_eigenvectorsOk = true;
557 }
558
559 m_nbrIterations = iparam[2];
560 m_nbrConverged = iparam[4];
561
562 m_info = Success;
563 }
564
565 delete[] select;
566 }
567
568 delete[] v;
569 delete[] iparam;
570 delete[] ipntr;
571 delete[] workd;
572 delete[] workl;
573 delete[] resid;
574
575 m_isInitialized = true;
576
577 return *this;
578}
579
580// Single precision
581//
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,
584 int *info);
585
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);
589
590// Double precision
591//
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,
594 int *info);
595
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);
599
600namespace internal {
601
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)
608 }
609
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,
613 int *ierr) {
614 EIGEN_STATIC_ASSERT(!NumTraits<Scalar>::IsComplex, NUMERIC_TYPE_MUST_BE_REAL)
615 }
616};
617
618template <>
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,
622 int *info) {
623 ssaupd_(ido, bmat, n, which, nev, tol, resid, ncv, v, ldv, iparam, ipntr, workd, workl, lworkl, info);
624 }
625
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);
631 }
632};
633
634template <>
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,
638 int *info) {
639 dsaupd_(ido, bmat, n, which, nev, tol, resid, ncv, v, ldv, iparam, ipntr, workd, workl, lworkl, info);
640 }
641
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);
647 }
648};
649
650template <typename MatrixSolver, typename MatrixType, typename Scalar, bool BisSPD>
651struct OP {
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);
654};
655
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) {
659 // OP = L^{-1} A L^{-T} (B = LL^T)
660 //
661 // First solve L^T out = in
662 //
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);
665
666 // Then compute out = A out
667 //
668 Matrix<Scalar, Dynamic, 1>::Map(out, n) = A * Matrix<Scalar, Dynamic, 1>::Map(out, n);
669
670 // Then solve L out = out
671 //
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));
674 }
675
676 static inline void project(MatrixSolver &OP, int n, int k, Scalar *vecs) {
677 // Solve L^T out = in
678 //
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);
683 }
684};
685
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...");
690 }
691
692 static inline void project(MatrixSolver &OP, int n, int k, Scalar *vecs) {
693 eigen_assert(false && "Should never be in here...");
694 }
695};
696
697} // end namespace internal
698
699} // end namespace Eigen
700
701#endif // EIGEN_ARPACKSELFADJOINTEIGENSOLVER_H
ComputationInfo
NumericalIssue
ComputeEigenvectors
Namespace containing all symbols from the Eigen library.
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index