1 #ifndef SYM_EIGS_SOLVER_H
2 #define SYM_EIGS_SOLVER_H
5 #include <Eigen/Eigenvalues>
14 #include "UpperHessenbergQR.h"
15 #include "MatOp/DenseGenMatProd.h"
16 #include "MatOp/DenseSymShiftSolve.h"
135 template <
typename Scalar = double,
141 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> Matrix;
142 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, 1> Vector;
143 typedef Eigen::Array<Scalar, Eigen::Dynamic, 1> Array;
144 typedef Eigen::Array<bool, Eigen::Dynamic, 1> BoolArray;
145 typedef Eigen::Map<const Matrix> MapMat;
146 typedef Eigen::Map<const Vector> MapVec;
147 typedef Eigen::SelfAdjointEigenSolver<Matrix> EigenSolver;
148 typedef std::pair<Scalar, int> SortPair;
182 void factorize_from(
int from_k,
int to_m,
const Vector &fk)
184 if(to_m <= from_k)
return;
188 Vector v(dim_n), w(dim_n);
189 Scalar beta = 0.0, Hii = 0.0;
191 fac_H.rightCols(ncv - from_k).setZero();
192 fac_H.block(from_k, 0, ncv - from_k, from_k).setZero();
193 for(
int i = from_k; i <= to_m - 1; i++)
196 v.noalias() = fac_f / beta;
198 fac_H(i, i - 1) = beta;
200 op->perform_op(v.data(), w.data());
204 fac_H(i - 1, i) = beta;
207 fac_f.noalias() = w - beta * fac_V.col(i - 1) - Hii * v;
212 Scalar v1f = fac_f.dot(fac_V.col(0));
213 if(v1f > prec || v1f < -prec)
216 Vf.tail(i) = fac_V.block(0, 1, dim_n, i).transpose() * fac_f;
218 fac_f -= fac_V.leftCols(i + 1) * Vf;
234 for(
int i = k; i < ncv; i++)
237 fac_H.diagonal().array() -= ritz_val[i];
246 fac_H.diagonal().array() += ritz_val[i];
251 Vector fk = fac_f * em[k - 1] + fac_V.col(k) * fac_H(k, k - 1);
252 factorize_from(k, ncv, fk);
257 int num_converged(Scalar tol)
260 Array thresh = tol * ritz_val.head(nev).array().abs().max(prec);
261 Array resid = ritz_vec.template bottomRows<1>().transpose().array().abs() * fac_f.norm();
263 ritz_conv = (resid < thresh);
265 return ritz_conv.cast<
int>().sum();
269 int nev_adjusted(
int nconv)
274 nev_new = nev + std::min(nconv, (ncv - nev) / 2);
275 if(nev == 1 && ncv >= 6)
277 else if(nev == 1 && ncv > 2)
284 void retrieve_ritzpair()
286 EigenSolver eig(fac_H);
287 Vector evals = eig.eigenvalues();
288 Matrix evecs = eig.eigenvectors();
290 std::vector<SortPair> pairs(ncv);
291 EigenvalueComparator<Scalar, SelectionRule> comp;
292 for(
int i = 0; i < ncv; i++)
294 pairs[i].first = evals[i];
297 std::sort(pairs.begin(), pairs.end(), comp);
308 std::vector<SortPair> pairs_copy(pairs);
309 for(
int i = 0; i < ncv; i++)
314 pairs[i] = pairs_copy[i / 2];
316 pairs[i] = pairs_copy[ncv - 1 - i / 2];
321 for(
int i = 0; i < ncv; i++)
323 ritz_val[i] = pairs[i].first;
325 for(
int i = 0; i < nev; i++)
327 ritz_vec.col(i) = evecs.col(pairs[i].second);
334 virtual void sort_ritzpair()
336 std::vector<SortPair> pairs(nev);
337 EigenvalueComparator<Scalar, LARGEST_MAGN> comp;
338 for(
int i = 0; i < nev; i++)
340 pairs[i].first = ritz_val[i];
343 std::sort(pairs.begin(), pairs.end(), comp);
345 Matrix new_ritz_vec(ncv, nev);
346 BoolArray new_ritz_conv(nev);
348 for(
int i = 0; i < nev; i++)
350 ritz_val[i] = pairs[i].first;
351 new_ritz_vec.col(i) = ritz_vec.col(pairs[i].second);
352 new_ritz_conv[i] = ritz_conv[pairs[i].second];
355 ritz_vec.swap(new_ritz_vec);
356 ritz_conv.swap(new_ritz_conv);
381 ncv(ncv_ > dim_n ? dim_n : ncv_),
384 prec(std::pow(std::numeric_limits<Scalar>::epsilon(), Scalar(2.0 / 3)))
386 if(nev_ < 1 || nev_ > dim_n - 1)
387 throw std::invalid_argument(
"nev must satisfy 1 <= nev <= n - 1, n is the size of matrix");
389 if(ncv_ <= nev_ || ncv_ > dim_n)
390 throw std::invalid_argument(
"ncv must satisfy nev < ncv <= n, n is the size of matrix");
402 void init(
const Scalar *init_resid)
405 fac_V.resize(dim_n, ncv);
406 fac_H.resize(ncv, ncv);
408 ritz_val.resize(ncv);
409 ritz_vec.resize(ncv, nev);
410 ritz_conv.resize(nev);
423 std::copy(init_resid, init_resid + dim_n, v.data());
424 Scalar vnorm = v.norm();
426 throw std::invalid_argument(
"initial residual vector cannot be zero");
430 op->perform_op(v.data(), w.data());
433 fac_H(0, 0) = v.dot(w);
434 fac_f = w - v * fac_H(0, 0);
447 Vector init_resid = Vector::Random(dim_n);
448 init_resid.array() -= 0.5;
449 init(init_resid.data());
460 int compute(
int maxit = 1000, Scalar tol = 1e-10)
463 factorize_from(1, ncv, fac_f);
466 int i, nconv, nev_adj;
467 for(i = 0; i < maxit; i++)
469 nconv = num_converged(tol);
473 nev_adj = nev_adjusted(nconv);
481 return std::min(nev, nconv);
503 int nconv = ritz_conv.cast<
int>().sum();
510 for(
int i = 0; i < nev; i++)
514 res[j] = ritz_val[i];
531 int nconv = ritz_conv.cast<
int>().sum();
532 Matrix res(dim_n, nconv);
537 Matrix ritz_vec_conv(ncv, nconv);
539 for(
int i = 0; i < nev; i++)
543 ritz_vec_conv.col(j) = ritz_vec.col(i);
548 res.noalias() = fac_V * ritz_vec_conv;
683 template <
typename Scalar = double,
689 typedef Eigen::Array<Scalar, Eigen::Dynamic, 1> Array;
696 Array ritz_val_org = Scalar(1.0) / this->ritz_val.head(this->nev).array() + sigma;
697 this->ritz_val.head(this->nev) = ritz_val_org;
720 SymEigsSolver<Scalar, SelectionRule, OpType>(op_, nev_, ncv_),
723 this->op->set_shift(sigma);
727 #endif // SYM_EIGS_SOLVER_H
void compute(ConstGenericMatrix &mat)
SymEigsSolver(OpType *op_, int nev_, int ncv_)
SymEigsShiftSolver(OpType *op_, int nev_, int ncv_, Scalar sigma_)
int compute(int maxit=1000, Scalar tol=1e-10)
void apply_YQ(GenericMatrix Y)
void init(const Scalar *init_resid)
void apply_QtY(Vector &Y)