votca 2024.2-dev
Loading...
Searching...
No Matches
davidsonsolver.h
Go to the documentation of this file.
1/*
2 * Copyright 2009-2020 The VOTCA Development Team (http://www.votca.org)
3 *
4 * Licensed under the Apache License, Version 2.0 (the "License");
5 * you may not use this file except in compliance with the License.
6 * You may obtain a copy of the License at
7 *
8 * http://www.apache.org/licenses/LICENSE-2.0
9 *
10 * Unless required by applicable law or agreed to in writing, software
11 * distributed under the License is distributed on an "AS IS" BASIS,
12 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
13 * See the License for the specific language governing permissions and
14 * limitations under the License.
15 *
16 */
17
18#pragma once
19#ifndef VOTCA_XTP_DAVIDSONSOLVER_H
20#define VOTCA_XTP_DAVIDSONSOLVER_H
21
22// Standard includes
23#include <chrono>
24#include <iostream>
25#include <stdexcept>
26
27// Third party includes
28#include <boost/format.hpp>
29
30// Local VOTCA includes
31#include "eigen.h"
32#include "logger.h"
33
34namespace votca {
35namespace xtp {
36
48
49 public:
51
52 void set_iter_max(Index N) { this->iter_max_ = N; }
54 void set_tolerance(std::string tol);
55 void set_correction(std::string method);
56 void set_size_update(std::string update_size);
57 void set_matrix_type(std::string mt);
58
59 Eigen::ComputationInfo info() const { return info_; }
60 Eigen::VectorXd eigenvalues() const { return this->eigenvalues_; }
61 Eigen::MatrixXd eigenvectors() const { return this->eigenvectors_; }
62 Index num_iterations() const { return this->i_iter_; }
63
64 template <typename MatrixReplacement>
65 void solve(const MatrixReplacement &A, Index neigen,
66 Index size_initial_guess = 0) {
67
68 if (max_search_space_ < neigen) {
69 max_search_space_ = neigen * 5;
70 }
71 std::chrono::time_point<std::chrono::system_clock> start =
72 std::chrono::system_clock::now();
73 Index op_size = A.rows();
74
75 checkOptions(op_size);
76 printOptions(op_size);
77
78 // initial guess size
79 if (size_initial_guess == 0) {
80 size_initial_guess = 2 * neigen;
81 }
82
83 restart_size_ = size_initial_guess;
84
85 // get the diagonal of the operator
86 this->Adiag_ = A.diagonal();
87
88 // target the lowest diagonal element
89 ProjectedSpace proj = initProjectedSpace(neigen, size_initial_guess);
90 RitzEigenPair rep;
92 << TimeStamp() << " iter\tSearch Space\tNorm" << std::flush;
93
94 for (i_iter_ = 0; i_iter_ < iter_max_; i_iter_++) {
95
96 updateProjection(A, proj);
97
98 rep = getRitzEigenPairs(proj);
99
100 bool converged = checkConvergence(rep, proj, neigen);
101
102 printIterationData(rep, proj, neigen);
103
104 bool last_iter = i_iter_ == (iter_max_ - 1);
105
106 if (converged) {
107 storeConvergedData(rep, neigen);
108 break;
109 } else if (last_iter) {
110 storeNotConvergedData(rep, proj.root_converged, neigen);
111 break;
112 }
113 Index extension_size = extendProjection(rep, proj);
114 bool do_restart = (proj.search_space() > max_search_space_);
115
116 if (do_restart) {
117 restart(rep, proj, extension_size);
118 }
119 }
120
121 printTiming(start);
122 }
123
124 private:
125 using ArrayXb = Eigen::Array<bool, Eigen::Dynamic, 1>;
129 double tol_ = 1E-4;
131 Eigen::VectorXd Adiag_;
133 enum CORR { DPR, OLSEN };
135
136 enum UPDATE { MIN, SAFE, MAX };
138
141
142 Eigen::VectorXd eigenvalues_;
143 Eigen::MatrixXd eigenvectors_;
144 Eigen::ComputationInfo info_ = Eigen::ComputationInfo::NoConvergence;
145
147 Eigen::VectorXd lambda; // eigenvalues
148 Eigen::MatrixXd q; // Ritz (or harmonic Ritz) eigenvectors
149 Eigen::MatrixXd U; // eigenvectors of the small subspace
150 Eigen::MatrixXd res; // residues of the pairs
151 Eigen::ArrayXd res_norm() const {
152 return res.colwise().norm();
153 } // norm of the residues
154 };
155
157 Eigen::MatrixXd V; // basis of vectors
158 Eigen::MatrixXd AV; // A * V
159 Eigen::MatrixXd T; // V.T * A * V
161 return V.cols();
162 }; // size of the projection i.e. number of cols in V
163 Index size_update; // size update ...
164 ArrayXb root_converged; // keep track of which root have onverged
165
166 // These are only used for harmonic ritz in the non-hermitian case
167 Eigen::MatrixXd AAV; // A*A*V
168 Eigen::MatrixXd B; // V.T *A*A*V
169 };
170
171 template <typename MatrixReplacement>
172 void updateProjection(const MatrixReplacement &A,
173 ProjectedSpace &proj) const {
174
175 if (i_iter_ == 0) {
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;
181 }
182
183 } else {
184 /* if we use a Gram Schmid(GS) orthogonalisation we do not have to
185 recompute the entire projection as GS doesn't change the original
186 subspace*/
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);
192
193 proj.T.conservativeResize(new_dim, new_dim);
194 proj.T.rightCols(nvec) = proj.V.transpose() * proj.AV.rightCols(nvec);
195
197 proj.T.bottomLeftCorner(nvec, old_dim) =
198 proj.T.topRightCorner(old_dim, nvec).transpose();
199
200 } else {
201 proj.T.bottomLeftCorner(nvec, old_dim) =
202 proj.V.rightCols(nvec).transpose() * proj.AV.leftCols(old_dim);
203
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);
210 }
211 }
212 }
213 RitzEigenPair getRitzEigenPairs(const ProjectedSpace &proj) const;
214
215 Eigen::MatrixXd qr(const Eigen::MatrixXd &A) const;
216 RitzEigenPair getHarmonicRitz(const ProjectedSpace &proj) const;
217
218 RitzEigenPair getRitz(const ProjectedSpace &proj) const;
219
220 Index getSizeUpdate(Index neigen) const;
221
222 void checkOptions(Index operator_size);
223
224 void printOptions(Index operator_size) const;
225
226 void printTiming(
227 const std::chrono::time_point<std::chrono::system_clock> &start) const;
228
229 void printIterationData(const RitzEigenPair &rep, const ProjectedSpace &proj,
230 Index neigen) const;
231
232 ArrayXl argsort(const Eigen::VectorXd &V) const;
233
234 Eigen::MatrixXd setupInitialEigenvectors(Index size_initial_guess) const;
235
236 Eigen::MatrixXd extract_vectors(const Eigen::MatrixXd &V,
237 const ArrayXl &idx) const;
238
239 void orthogonalize(Eigen::MatrixXd &V, Index nupdate) const;
240 void gramschmidt(Eigen::MatrixXd &A, Index nstart) const;
241
242 Eigen::VectorXd computeCorrectionVector(const Eigen::VectorXd &qj,
243 double lambdaj,
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;
248
249 ProjectedSpace initProjectedSpace(Index neigen,
250 Index size_initial_guess) const;
251
252 Index extendProjection(const RitzEigenPair &rep, ProjectedSpace &proj) const;
253
254 bool checkConvergence(const RitzEigenPair &rep, ProjectedSpace &proj,
255 Index neigen) const;
256
257 void restart(const RitzEigenPair &rep, ProjectedSpace &proj,
258 Index newtestvectors) const;
259
260 void storeConvergedData(const RitzEigenPair &rep, Index neigen);
261
262 void storeNotConvergedData(const RitzEigenPair &rep,
263 const ArrayXb &root_converged, Index neigen);
264
265 void storeEigenPairs(const RitzEigenPair &rep, Index neigen);
266};
267
268} // namespace xtp
269} // namespace votca
270
271#endif // VOTCA_XTP_DAVIDSONSOLVER_H
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
ArrayXl argsort(const Eigen::VectorXd &V) const
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 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
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::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.
Definition logger.h:164
Timestamp returns the current time as a string Example: cout << TimeStamp()
Definition logger.h:224
#define XTP_LOG(level, log)
Definition logger.h:40
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26
Eigen::Array< votca::Index, Eigen::Dynamic, 1 > ArrayXl
Definition eigen.h:36