19#ifndef VOTCA_XTP_DAVIDSONSOLVER_H
20#define VOTCA_XTP_DAVIDSONSOLVER_H
28#include <boost/format.hpp>
59 Eigen::ComputationInfo
info()
const {
return info_; }
64 template <
typename MatrixReplacement>
66 Index size_initial_guess = 0) {
71 std::chrono::time_point<std::chrono::system_clock> start =
72 std::chrono::system_clock::now();
73 Index op_size = A.rows();
79 if (size_initial_guess == 0) {
80 size_initial_guess = 2 * neigen;
86 this->
Adiag_ = A.diagonal();
92 <<
TimeStamp() <<
" iter\tSearch Space\tNorm" << std::flush;
109 }
else if (last_iter) {
117 restart(rep, proj, extension_size);
125 using ArrayXb = Eigen::Array<bool, Eigen::Dynamic, 1>;
144 Eigen::ComputationInfo
info_ = Eigen::ComputationInfo::NoConvergence;
152 return res.colwise().norm();
171 template <
typename MatrixReplacement>
176 proj.
AV = A * proj.
V;
177 proj.
T = proj.
V.transpose() * proj.
AV;
179 proj.
AAV = A * proj.
AV;
180 proj.
B = proj.
V.transpose() * proj.
AAV;
187 Index old_dim = proj.
AV.cols();
188 Index new_dim = proj.
V.cols();
189 Index nvec = new_dim - old_dim;
190 proj.
AV.conservativeResize(Eigen::NoChange, new_dim);
191 proj.
AV.rightCols(nvec) = A * proj.
V.rightCols(nvec);
193 proj.
T.conservativeResize(new_dim, new_dim);
194 proj.
T.rightCols(nvec) = proj.
V.transpose() * proj.
AV.rightCols(nvec);
197 proj.
T.bottomLeftCorner(nvec, old_dim) =
198 proj.
T.topRightCorner(old_dim, nvec).transpose();
201 proj.
T.bottomLeftCorner(nvec, old_dim) =
202 proj.
V.rightCols(nvec).transpose() * proj.
AV.leftCols(old_dim);
204 proj.
AAV.conservativeResize(Eigen::NoChange, new_dim);
205 proj.
AAV.rightCols(nvec) = A * proj.
AV.rightCols(nvec);
206 proj.
B.conservativeResize(new_dim, new_dim);
207 proj.
B.rightCols(nvec) = proj.
V.transpose() * proj.
AAV.rightCols(nvec);
208 proj.
B.bottomLeftCorner(nvec, old_dim) =
209 proj.
V.rightCols(nvec).transpose() * proj.
AAV.leftCols(old_dim);
215 Eigen::MatrixXd
qr(
const Eigen::MatrixXd &A)
const;
218 RitzEigenPair
getRitz(
const ProjectedSpace &proj)
const;
227 const std::chrono::time_point<std::chrono::system_clock> &start)
const;
244 const Eigen::VectorXd &Aqj)
const;
245 Eigen::VectorXd
dpr(
const Eigen::VectorXd &r,
double lambda)
const;
246 Eigen::VectorXd
olsen(
const Eigen::VectorXd &r,
const Eigen::VectorXd &x,
247 double lambda)
const;
250 Index size_initial_guess)
const;
257 void restart(
const RitzEigenPair &rep, ProjectedSpace &proj,
258 Index newtestvectors)
const;
Use Davidson algorithm to solve A*V=E*V.
Index getSizeUpdate(Index neigen) const
RitzEigenPair getRitzEigenPairs(const ProjectedSpace &proj) const
bool checkConvergence(const RitzEigenPair &rep, ProjectedSpace &proj, Index neigen) const
void restart(const RitzEigenPair &rep, ProjectedSpace &proj, Index newtestvectors) const
ProjectedSpace initProjectedSpace(Index neigen, Index size_initial_guess) const
void printOptions(Index operator_size) const
Eigen::MatrixXd eigenvectors_
ArrayXl argsort(const Eigen::VectorXd &V) const
DavidsonSolver(Logger &log)
Index extendProjection(const RitzEigenPair &rep, ProjectedSpace &proj) const
void set_max_search_space(Index N)
Eigen::VectorXd eigenvalues() const
void solve(const MatrixReplacement &A, Index neigen, Index size_initial_guess=0)
void set_iter_max(Index N)
void checkOptions(Index operator_size)
void set_matrix_type(std::string mt)
Eigen::MatrixXd setupInitialEigenvectors(Index size_initial_guess) const
void set_correction(std::string method)
Eigen::ComputationInfo info() const
void set_size_update(std::string update_size)
RitzEigenPair getHarmonicRitz(const ProjectedSpace &proj) const
Index num_iterations() const
void set_tolerance(std::string tol)
Eigen::MatrixXd eigenvectors() const
Eigen::VectorXd dpr(const Eigen::VectorXd &r, double lambda) const
void gramschmidt(Eigen::MatrixXd &A, Index nstart) const
Eigen::VectorXd computeCorrectionVector(const Eigen::VectorXd &qj, double lambdaj, const Eigen::VectorXd &Aqj) const
Eigen::MatrixXd extract_vectors(const Eigen::MatrixXd &V, const ArrayXl &idx) const
void printTiming(const std::chrono::time_point< std::chrono::system_clock > &start) const
void storeEigenPairs(const RitzEigenPair &rep, Index neigen)
Eigen::Array< bool, Eigen::Dynamic, 1 > ArrayXb
RitzEigenPair getRitz(const ProjectedSpace &proj) const
void storeConvergedData(const RitzEigenPair &rep, Index neigen)
void orthogonalize(Eigen::MatrixXd &V, Index nupdate) const
void updateProjection(const MatrixReplacement &A, ProjectedSpace &proj) const
Eigen::ComputationInfo info_
void printIterationData(const RitzEigenPair &rep, const ProjectedSpace &proj, Index neigen) const
Eigen::VectorXd olsen(const Eigen::VectorXd &r, const Eigen::VectorXd &x, double lambda) const
Eigen::VectorXd eigenvalues_
CORR davidson_correction_
Eigen::MatrixXd qr(const Eigen::MatrixXd &A) const
void storeNotConvergedData(const RitzEigenPair &rep, const ArrayXb &root_converged, Index neigen)
Logger is used for thread-safe output of messages.
Timestamp returns the current time as a string Example: cout << TimeStamp()
#define XTP_LOG(level, log)
base class for all analysis tools
Index search_space() const
Eigen::ArrayXd res_norm() const
Eigen::Array< votca::Index, Eigen::Dynamic, 1 > ArrayXl