lambda-lanczos 2.1.1
Loading...
Searching...
No Matches
lambda_lanczos.hpp
Go to the documentation of this file.
1#ifndef LAMBDA_LANCZOS_H_
2#define LAMBDA_LANCZOS_H_
3
4#include <algorithm>
5#include <cassert>
6#include <cmath>
7#include <functional>
8#include <iostream>
9#include <limits>
10#include <numeric>
11#include <random>
12#include <tuple>
13#include <utility>
14#include <vector>
15
16#include "eigenpair_manager.hpp"
19
20namespace lambda_lanczos {
32template <typename T, typename LT>
33inline auto compute_eigenvectors(const std::vector<T>& alpha,
34 const std::vector<T>& beta,
35 const std::vector<std::vector<LT>>& u,
36 const bool find_maximum,
37 const size_t num_of_eigenvalues) -> std::vector<std::vector<decltype(T() + LT())>> {
38 const auto m = alpha.size();
39 const auto n = u[0].size();
40
41 std::vector<T> tridiagonal_eigenvalues;
42 std::vector<std::vector<T>> tridiagonal_eigenvectors;
43
44 tridiagonal::tridiagonal_eigenpairs(alpha, beta, tridiagonal_eigenvalues, tridiagonal_eigenvectors);
45
46 std::vector<std::vector<LT>> eigenvectors;
47 for (size_t i = 0; i < num_of_eigenvalues; ++i) {
48 eigenvectors.emplace_back(n);
49 }
50
51 for (size_t index = 0; index < num_of_eigenvalues; index++) {
52 size_t index_tri = find_maximum ? m - index - 1 : index;
53 for (size_t k = m; k-- > 0;) {
54 for (size_t i = 0; i < n; ++i) {
55 eigenvectors[index][i] += tridiagonal_eigenvectors[index_tri][k] * u[k][i];
56 }
57 }
58 util::normalize(eigenvectors[index]);
59 }
60
61 return eigenvectors;
62}
63
70template <typename T>
72 public:
79 static void init(std::vector<T>& v) {
80 std::random_device dev;
81 std::mt19937 mt(dev());
82 std::uniform_real_distribution<T> rand((T)(-1.0), (T)(1.0));
83
84 size_t n = v.size();
85 for (size_t i = 0; i < n; ++i) {
86 v[i] = rand(mt);
87 }
88 }
89};
90
91template <typename T>
92struct VectorRandomInitializer<std::complex<T>> {
93 public:
94 static void init(std::vector<std::complex<T>>& v) {
95 std::random_device dev;
96 std::mt19937 mt(dev());
97 std::uniform_real_distribution<T> rand((T)(-1.0), (T)(1.0));
98
99 size_t n = v.size();
100 for (size_t i = 0; i < n; ++i) {
101 v[i] = std::complex<T>(rand(mt), rand(mt));
102 }
103 }
104};
105
109template <typename T>
111 private:
115 template <typename n_type>
117
118 public:
126 std::function<void(const std::vector<T>& in, std::vector<T>& out)> mv_mul;
127
133 std::function<void(std::vector<T>& vec)> init_vector = VectorRandomInitializer<T>::init;
134
150 real_t<T> eps = std::numeric_limits<real_t<T>>::epsilon() * 1e3;
151
154
156 size_t num_eigs = 1;
157
166
174
182
183 private:
187 std::vector<size_t> iter_counts;
188
189 public:
200 LambdaLanczos(std::function<void(const std::vector<T>&, std::vector<T>&)> mv_mul,
201 size_t matrix_size,
202 bool find_maximum,
203 size_t num_eigs)
204 : mv_mul(mv_mul),
209
216 template <typename Iterable>
217 size_t run_iteration(std::vector<real_t<T>>& eigvalues,
218 std::vector<std::vector<T>>& eigvecs,
219 size_t nroot,
220 Iterable orthogonalizeTo) const {
221 std::vector<std::vector<T>> u; // Lanczos vectors
222 std::vector<real_t<T>> alpha; // Diagonal elements of an approximated tridiagonal matrix
223 std::vector<real_t<T>> beta; // Subdiagonal elements of an approximated tridiagonal matrix
224
225 u.reserve(this->initial_vector_size);
226 alpha.reserve(this->initial_vector_size);
227 beta.reserve(this->initial_vector_size);
228
229 const auto n = this->matrix_size;
230
231 u.emplace_back(n);
232 this->init_vector(u[0]);
233 util::schmidt_orth(u[0], orthogonalizeTo.cbegin(), orthogonalizeTo.cend());
234 util::normalize(u[0]);
235
236 std::vector<real_t<T>> evs; // Calculated eigenvalue
237 std::vector<real_t<T>> pevs; // Previous eigenvalue
238
239 size_t itern = this->max_iteration;
240 for (size_t k = 1; k <= this->max_iteration; ++k) {
241 /* au = (A + offset*E)uk, here E is the identity matrix */
242 std::vector<T> au(n, 0.0); // Temporal storage to store matrix-vector multiplication result
243 this->mv_mul(u[k - 1], au);
244 for (size_t i = 0; i < n; ++i) {
245 au[i] += u[k - 1][i] * this->eigenvalue_offset;
246 }
247
248 alpha.push_back(std::real(util::inner_prod(u[k - 1], au)));
249
250 u.push_back(std::move(au));
251 for (size_t i = 0; i < n; ++i) {
252 if (k == 1) {
253 u[k][i] = u[k][i] - alpha[k - 1] * u[k - 1][i];
254 } else {
255 u[k][i] = u[k][i] - beta[k - 2] * u[k - 2][i] - alpha[k - 1] * u[k - 1][i];
256 }
257 }
258
259 util::schmidt_orth(u[k], orthogonalizeTo.cbegin(), orthogonalizeTo.cend());
260 util::schmidt_orth(u[k], u.begin(), u.end() - 1);
261
262 beta.push_back(util::norm(u[k]));
263
264 size_t num_eigs_to_calculate = std::min(nroot, alpha.size());
265 evs = std::vector<real_t<T>>();
266
267 std::vector<real_t<T>> eigvals_all(alpha.size());
268 tridiagonal::tridiagonal_eigenvalues(alpha, beta, eigvals_all);
269 if (this->find_maximum) {
270 for (size_t i = 0; i < num_eigs_to_calculate; ++i) {
271 evs.push_back(eigvals_all[eigvals_all.size() - i - 1]);
272 }
273 } else {
274 for (size_t i = 0; i < num_eigs_to_calculate; ++i) {
275 evs.push_back(eigvals_all[i]);
276 }
277 }
278
279 const real_t<T> zero_threshold = std::numeric_limits<real_t<T>>::epsilon() * 1e1;
280 if (beta.back() < zero_threshold) {
281 itern = k;
282 break;
283 }
284
285 util::normalize(u[k]);
286
287 /*
288 * Only break loop if convergence condition is met for all requied roots
289 */
290 bool break_cond = true;
291 if (pevs.size() != evs.size()) {
292 break_cond = false;
293 } else {
294 for (size_t iroot = 0; iroot < nroot; ++iroot) {
295 const auto& ev = evs[iroot];
296 const auto& pev = pevs[iroot];
297 if (std::abs(ev - pev) >= std::min(std::abs(ev), std::abs(pev)) * this->eps) {
298 break_cond = false;
299 break;
300 }
301 }
302 }
303
304 if (break_cond) {
305 itern = k;
306 break;
307 } else {
308 pevs = evs;
309 }
310 }
311
312 eigvalues = evs;
313 eigvecs.resize(eigvalues.size());
314 beta.back() = 0.0;
315
316 eigvecs = compute_eigenvectors(alpha, beta, u, find_maximum, eigvalues.size());
317 for (size_t i = 0; i < eigvalues.size(); ++i) {
318 eigvalues[i] -= this->eigenvalue_offset;
319 }
320
321 return itern;
322 }
323
330 void run(std::vector<real_t<T>>& eigenvalues, std::vector<std::vector<T>>& eigenvectors) {
331 this->iter_counts = std::vector<size_t>();
333
334 while (true) {
335 std::vector<real_t<T>> eigenvalues_current;
336 std::vector<std::vector<T>> eigenvectors_current;
337
338 size_t nroot = std::min(num_eigs_per_iteration, this->matrix_size - ep_manager.size());
339
340 size_t iter_count =
341 this->run_iteration(eigenvalues_current, eigenvectors_current, nroot, ep_manager.getEigenvectors());
342 this->iter_counts.push_back(iter_count);
343
344 bool nothing_added = ep_manager.insertEigenpairs(eigenvalues_current, eigenvectors_current);
345
346 if (nothing_added) {
347 break;
348 }
349
350 if (this->num_eigs == 1) {
351 // If only one eigenpair is required, iteration can be stopped
352 break;
353 }
354 }
355
356 const auto& eigenpairs = ep_manager.getEigenpairs();
357 eigenvalues = std::vector<real_t<T>>();
358 eigenvalues.reserve(eigenpairs.size());
359 eigenvectors = std::vector<std::vector<T>>();
360 eigenvectors.reserve(eigenvectors.size());
361
362 for (auto& eigenpair : eigenpairs) {
363 eigenvalues.push_back(std::move(eigenpair.first));
364 eigenvectors.push_back(std::move(eigenpair.second));
365 }
366 }
367
376 std::tuple<std::vector<real_t<T>>, std::vector<std::vector<T>>> run() {
377 std::vector<real_t<T>> eigenvalues;
378 std::vector<std::vector<T>> eigenvectors;
379 for (size_t i = 0; i < this->num_eigs; ++i) {
380 eigenvectors.emplace_back(this->matrix_size);
381 }
382
383 this->run(eigenvalues, eigenvectors);
384
385 return {eigenvalues, eigenvectors};
386 }
387
394 void run(real_t<T>& eigenvalue, std::vector<T>& eigenvector) {
395 const size_t num_eigs_tmp = this->num_eigs;
396 this->num_eigs = 1;
397
398 std::vector<real_t<T>> eigenvalues(1);
399 std::vector<std::vector<T>> eigenvectors(1);
400
401 this->run(eigenvalues, eigenvectors);
402
403 this->num_eigs = num_eigs_tmp;
404
405 eigenvalue = eigenvalues[0];
406 eigenvector = std::move(eigenvectors[0]);
407 }
408
412 const std::vector<size_t>& getIterationCounts() const {
413 return iter_counts;
414 }
415};
416
417} /* namespace lambda_lanczos */
418
419#endif /* LAMBDA_LANCZOS_H_ */
size_t max_iteration
Iteration limit of Lanczos algorithm, set to matrix_size automatically.
Definition lambda_lanczos.hpp:138
std::tuple< std::vector< real_t< T > >, std::vector< std::vector< T > > > run()
Executes Lanczos algorithm and return result as a tuple.
Definition lambda_lanczos.hpp:376
size_t num_eigs_per_iteration
(Not necessary to change)
Definition lambda_lanczos.hpp:173
util::real_t< n_type > real_t
See util::real_t for details.
Definition lambda_lanczos.hpp:116
real_t< T > eps
Convergence threshold of Lanczos iteration.
Definition lambda_lanczos.hpp:150
LambdaLanczos(std::function< void(const std::vector< T > &, std::vector< T > &)> mv_mul, size_t matrix_size, bool find_maximum, size_t num_eigs)
Constructs Lanczos calculation engine.
Definition lambda_lanczos.hpp:200
std::function< void(const std::vector< T > &in, std::vector< T > &out)> mv_mul
Matrix-vector multiplication routine.
Definition lambda_lanczos.hpp:126
void run(std::vector< real_t< T > > &eigenvalues, std::vector< std::vector< T > > &eigenvectors)
Executes Lanczos algorithm and stores the result into reference variables passed as arguments.
Definition lambda_lanczos.hpp:330
size_t run_iteration(std::vector< real_t< T > > &eigvalues, std::vector< std::vector< T > > &eigvecs, size_t nroot, Iterable orthogonalizeTo) const
Not documented (In most cases, run() is preferred).
Definition lambda_lanczos.hpp:217
std::function< void(std::vector< T > &vec)> init_vector
Function to initialize the initial Lanczos vector.
Definition lambda_lanczos.hpp:133
size_t num_eigs
Number of eigenpairs to be calculated.
Definition lambda_lanczos.hpp:156
std::vector< size_t > iter_counts
Iteration counts of the latest run.
Definition lambda_lanczos.hpp:187
size_t matrix_size
Dimension of the matrix to be diagonalized.
Definition lambda_lanczos.hpp:136
real_t< T > eigenvalue_offset
Shifts the eigenvalues of the given matrix A.
Definition lambda_lanczos.hpp:165
const std::vector< size_t > & getIterationCounts() const
Returns the latest iteration counts.
Definition lambda_lanczos.hpp:412
bool find_maximum
true to calculate maximum eigenvalue, false to calculate minimum one.
Definition lambda_lanczos.hpp:153
size_t initial_vector_size
(Not necessary to change)
Definition lambda_lanczos.hpp:181
void run(real_t< T > &eigenvalue, std::vector< T > &eigenvector)
Executes Lanczos algorithm that calculate one eigenpair regardless of num_eigs.
Definition lambda_lanczos.hpp:394
Class to manage calculated eigenpairs.
Definition eigenpair_manager.hpp:22
size_t size() const
Definition eigenpair_manager.hpp:48
decltype(eigenpairs) & getEigenpairs()
Definition eigenpair_manager.hpp:77
lambda_lanczos::util::MapValueIterable< decltype(eigenpairs)> getEigenvectors() const
Definition eigenpair_manager.hpp:73
bool insertEigenpairs(std::vector< real_t< T > > &eigenvalues, std::vector< std::vector< T > > &eigenvectors)
Definition eigenpair_manager.hpp:52
size_t tridiagonal_eigenvalues(const std::vector< T > &alpha, const std::vector< T > &beta, std::vector< T > &eigenvalues)
Computes all eigenvalues for given tri-diagonal matrix using the Implicitly Shifted QR algorithm.
Definition lambda_lanczos_tridiagonal_impl.hpp:356
size_t tridiagonal_eigenpairs(const std::vector< T > &alpha, const std::vector< T > &beta, std::vector< T > &eigenvalues, std::vector< std::vector< T > > &eigenvectors, bool compute_eigenvector=true)
Computes all eigenpairs (eigenvalues and eigenvectors) for given tri-diagonal matrix using the Implic...
Definition lambda_lanczos_tridiagonal_impl.hpp:291
void normalize(std::vector< T > &vec)
Normalizes given vector.
Definition linear_algebra.hpp:78
real_t< T > norm(const std::vector< T > &vec)
Returns Euclidean norm of given vector.
Definition linear_algebra.hpp:57
typename realTypeMap< T >::type real_t
Type mapper from T to real type of T.
Definition common.hpp:102
T inner_prod(const std::vector< T > &v1, const std::vector< T > &v2)
Returns "mathematical" inner product of v1 and v2.
Definition linear_algebra.hpp:30
void schmidt_orth(std::vector< T > &uorth, ForwardIterator first, ForwardIterator last)
Orthogonalizes vector uorth with respect to orthonormal vectors defined by given iterators.
Definition linear_algebra.hpp:133
Definition eigenpair_manager.hpp:9
auto compute_eigenvectors(const std::vector< T > &alpha, const std::vector< T > &beta, const std::vector< std::vector< LT > > &u, const bool find_maximum, const size_t num_of_eigenvalues) -> std::vector< std::vector< decltype(T()+LT())> >
Computes the eigenvectors from Krylov subspace information.
Definition lambda_lanczos.hpp:33
static void init(std::vector< std::complex< T > > &v)
Definition lambda_lanczos.hpp:94
Template class to implement random vector initializer.
Definition lambda_lanczos.hpp:71
static void init(std::vector< T > &v)
Initialize given vector randomly in the range of [-1, 1].
Definition lambda_lanczos.hpp:79