lambda-lanczos 2.0.0
Loading...
Searching...
No Matches
exponentiator.hpp
Go to the documentation of this file.
1#ifndef LAMBDA_LANCZOS_EXPONENTIATOR_H_
2#define LAMBDA_LANCZOS_EXPONENTIATOR_H_
3
4#include <iostream>
5#include <iomanip>
6#include <vector>
7#include <tuple>
8#include <functional>
9#include <cassert>
10#include <limits>
11#include <cmath>
12#include <numeric>
13#include <random>
14#include <algorithm>
15#include "lambda_lanczos.hpp"
17
18
19namespace lambda_lanczos {
20
21
25template <typename T>
27private:
31 template <typename n_type>
33
34public:
42 std::function<void(const std::vector<T>& in, std::vector<T>& out)> mv_mul;
43
59 real_t<T> eps = std::numeric_limits<real_t<T>>::epsilon() * 1e2;
60
64 bool full_orthogonalize = false;
65
72 size_t initial_vector_size = 200;
73
81 Exponentiator(std::function<void(const std::vector<T>&, std::vector<T>&)> mv_mul, size_t matrix_size) :
83
88 size_t run(const T& a,
89 const std::vector<T>& input,
90 std::vector<T>& output) const {
91 assert(input.size() == this->matrix_size);
92
93 std::vector<std::vector<T>> u; // Lanczos vectors
94 std::vector<real_t<T>> alpha; // Diagonal elements of an approximated tridiagonal matrix
95 std::vector<real_t<T>> beta; // Subdiagonal elements of an approximated tridiagonal matrix
96
97 u.reserve(this->initial_vector_size);
98 alpha.reserve(this->initial_vector_size);
99 beta.reserve(this->initial_vector_size);
100
101 const auto n = this->matrix_size;
102
103 u.push_back(input);
104 util::normalize(u[0]);
105
106 std::vector<T> coeff_prev;
107
108 size_t itern = this->max_iteration;
109 for(size_t k = 1; k <= this->max_iteration; ++k) {
110 u.emplace_back(n, 0.0);
111 this->mv_mul(u[k-1], u[k]);
112
113 alpha.push_back(std::real(util::inner_prod(u[k-1], u[k])));
114
115 for(size_t i = 0; i < n; ++i) {
116 if(k == 1) {
117 u[k][i] = u[k][i] - alpha[k-1]*u[k-1][i];
118 } else {
119 u[k][i] = u[k][i] - beta[k-2]*u[k-2][i] - alpha[k-1]*u[k-1][i];
120 }
121 }
122
123 if(this->full_orthogonalize) {
124 util::schmidt_orth(u[k], u.begin(), u.end()-1);
125 }
126
127 std::vector<real_t<T>> ev(alpha.size());
128 std::vector<std::vector<real_t<T>>> p(alpha.size());
130
131 std::vector<T> coeff(alpha.size(), 0.0);
132 for(size_t i = 0; i < alpha.size(); ++i) {
133 for(size_t j = 0; j < alpha.size(); ++j) {
134 coeff[i] += p[j][i] * std::exp(a*ev[j]) * p[j][0];
135 }
136 }
137
138 /*std::cout << "beta:" << std::endl;
139 for(size_t i = 0; i < beta.size(); ++i) {
140 std::cout << i << " : " << beta[i] << std::endl;
141 }
142
143 std::cout << "coeff:" << std::endl;
144 for(size_t i = 0; i < coeff.size(); ++i) {
145 std::cout << i << " : " << coeff[i] << std::endl;
146 }*/
147
148 beta.push_back(util::norm(u[k]));
149
150
151 T overlap = 0.0;
152 for(size_t i = 0; i < coeff_prev.size(); ++i) {
153 overlap += util::typed_conj(coeff_prev[i])*coeff[i]; // Last element of coeff is not used here
154 }
155
156 coeff_prev = std::move(coeff);
157
158 const real_t<T> beta_threshold = std::numeric_limits<real_t<T>>::epsilon();
159 if(std::abs(1-std::abs(overlap)) < eps || beta.back() < beta_threshold) {
160 itern = k;
161 break;
162 }
163
164 util::normalize(u[k]);
165 }
166
167 output.resize(n);
168 std::fill(output.begin(), output.end(), T());
169 const T norm = util::norm(input);
170 for(size_t l = 0;l < coeff_prev.size(); ++l) {
171 for(size_t i = 0;i < n; ++i) {
172 output[i] += norm*coeff_prev[l]*u[l][i];
173 }
174 }
175
176 return itern;
177 }
178
179
180 size_t taylor_run(const T& a,
181 const std::vector<T>& input,
182 std::vector<T>& output) {
183 const size_t n = this->matrix_size;
184 assert(this->matrix_size == input.size());
185
186
187 if(a == T()) { // Zero check
188 output = input;
189 return 1;
190 }
191
192 std::vector<std::vector<T>> taylors;
193 taylors.push_back(input);
194
195 T factor = 1.0;
196 for(size_t k = 1; ; ++k) {
197 factor *= a/(T)k;
198 taylors.emplace_back(n, 0.0);
199 mv_mul(taylors[k-1], taylors[k]);
200
201 if(lambda_lanczos::util::norm(taylors[k]) * std::abs(factor) < eps) {
202 break;
203 }
204 }
205
206
207 /* Sum Taylor series backward */
208 output.resize(n);
209 std::fill(output.begin(), output.end(), 0.0);
210 for(size_t k = taylors.size(); k-- > 0; ) {
211 for(size_t i = 0; i < n; ++i) {
212 output[i] += taylors[k][i] * factor;
213 }
214
215 factor *= (T)k / a;
216 }
217
218 return taylors.size();
219 }
220};
221
222
223} /* namespace lambda_lanczos */
224
225#endif /* LAMBDA_LANCZOS_EXPONENTIATOR_H_ */
Calculation engine for Lanczos exponentiation.
Definition: exponentiator.hpp:26
size_t taylor_run(const T &a, const std::vector< T > &input, std::vector< T > &output)
Definition: exponentiator.hpp:180
real_t< T > eps
Convergence threshold of Lanczos iteration.
Definition: exponentiator.hpp:59
size_t max_iteration
Iteration limit of Lanczos algorithm, set to matrix_size automatically.
Definition: exponentiator.hpp:47
size_t matrix_size
Dimension of the matrix to be exponentiated.
Definition: exponentiator.hpp:45
std::function< void(const std::vector< T > &in, std::vector< T > &out)> mv_mul
Matrix-vector multiplication routine.
Definition: exponentiator.hpp:42
size_t initial_vector_size
(Not necessary to change)
Definition: exponentiator.hpp:72
util::real_t< n_type > real_t
See util::real_t for details.
Definition: exponentiator.hpp:32
size_t run(const T &a, const std::vector< T > &input, std::vector< T > &output) const
Apply matrix exponentiation exp(a*A) to input and store the result into output.
Definition: exponentiator.hpp:88
Exponentiator(std::function< void(const std::vector< T > &, std::vector< T > &)> mv_mul, size_t matrix_size)
Constructs Lanczos exponetiation engine.
Definition: exponentiator.hpp:81
bool full_orthogonalize
Flag to execute explicit Lanczos-vector orthogonalization.
Definition: exponentiator.hpp:64
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:302
void normalize(std::vector< T > &vec)
Normalizes given vector.
Definition: lambda_lanczos_util.hpp:180
real_t< T > norm(const std::vector< T > &vec)
Returns Euclidean norm of given vector.
Definition: lambda_lanczos_util.hpp:161
T typed_conj(const T &val)
Complex conjugate with type. This function returns the argument itself for real type,...
Definition: lambda_lanczos_util.hpp:135
T inner_prod(const std::vector< T > &v1, const std::vector< T > &v2)
Returns "mathematical" inner product of v1 and v2.
Definition: lambda_lanczos_util.hpp:148
void schmidt_orth(std::vector< T > &uorth, ForwardIterator first, ForwardIterator last)
Orthogonalizes vector uorth with respect to orthonormal vectors defined by given iterators.
Definition: lambda_lanczos_util.hpp:206
typename realTypeMap< T >::type real_t
Type mapper from T to real type of T.
Definition: lambda_lanczos_util.hpp:103
Definition: eigenpair_manager.hpp:10