-
Notifications
You must be signed in to change notification settings - Fork 2
/
krylov.hpp
341 lines (284 loc) · 9.95 KB
/
krylov.hpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
#ifndef __KRYLOV_HPP__
#define __KRYLOV_HPP__
#include <omp.h>
#include <Eigen/Core>
#include <Eigen/Dense>
#include <Eigen/Eigenvalues>
#include <Eigen/Sparse>
#include <fstream>
#include <iostream>
#include <random>
#include <unordered_map>
#include <unsupported/Eigen/MatrixFunctions>
/**
For convenience, we introduce some new types
*/
using basis = std::unordered_map<std::string, int>;
using inversebasis = std::unordered_map<int, std::string>;
using observable = Eigen::SparseMatrix<std::complex<double>, Eigen::RowMajor>;
using state = Eigen::VectorXcd;
using complex = std::complex<double>;
typedef Eigen::Triplet<complex> T;
/**
Converts an integer to a binary string representation, and
pads it until it has length L. This was done initially with
bitsets, but they can't be initialized with dynamically known
lengths.
*/
std::string int_to_bin(int i, int L) {
std::string binary_representation = "";
while (i > 0) {
binary_representation = std::to_string(i % 2) + binary_representation;
i /= 2;
}
while (binary_representation.length() < L)
binary_representation = "0" + binary_representation;
return binary_representation;
}
/**
Count the number of non-zero entries in the bitstring s
*/
int count_ones(std::string bitstring) {
int count = 0;
for (int i = 0; i < bitstring.length(); ++i)
if (bitstring[i] == '1') count++;
return count;
}
/**
Builds the allowed basis states for L sites and k particles
*/
basis build_basis(const int L, int k, inversebasis &inversebasis) {
int counter = 0;
basis basis;
for (int i = 0; i < pow(2, L); ++i) {
std::string binrep = int_to_bin(i, L);
if (count_ones(binrep) == k) {
inversebasis.insert({counter, binrep});
basis.insert({binrep, counter++});
}
}
return basis;
}
/**
Constructs the Hamiltonian matrix in sparse format
*/
observable build_hamiltonian(basis basis, inversebasis inversebasis, complex J,
complex F, complex U, complex W, int seed) {
int num_basis = (int)basis.size();
// Initialize the random number generator
std::default_random_engine generator(seed);
// And create a uniform distribution we can sample from
std::uniform_real_distribution<double> distribution(-1.0, 1.0);
std::vector<double> onsitedisorder;
int systemsize = inversebasis[0].length();
for (int i = 0; i < systemsize; ++i)
onsitedisorder.push_back(distribution(generator));
// Vector storing all the triplet pairs for setting the Sparse
// Hamiltonian later on.
std::vector<T> total;
// For each basis state..
#pragma omp parallel
{
// Each thread has its own local triplet list
std::vector<T> tripletList;
#pragma omp for nowait
for (int i = 0; i < num_basis; ++i) {
// Extract the string describing this state
std::string this_state = inversebasis[i];
//
// .. check for hopping
//
for (int site = 0; site < systemsize - 1; ++site) {
// Make a copy of the state
std::string newstate(this_state);
if (this_state[site] == '0' && this_state[site + 1] == '1') {
newstate[site] = '1';
newstate[site + 1] = '0';
basis::const_iterator got = basis.find(newstate);
int j = (int)got->second;
// Add forward and backward hopping entries
tripletList.push_back(T(i, j, J));
tripletList.push_back(T(j, i, std::conj(J)));
}
}
//
// .. then the local field and disorder
//
for (int site = 0; site < systemsize; ++site) {
complex localspin = (int)(this_state[site] - '0');
tripletList.push_back(T(i, i,
onsitedisorder.at(site) * localspin * W +
localspin * F * complex(site, 0)));
} // Local field and disorder
//
// .. and then the interactions
//
for (int site = 0; site < systemsize - 1; ++site) {
complex localspin = (int)(this_state[site] - '0');
complex nextspin = (int)(this_state[site + 1] - '0');
tripletList.push_back(T(i, i, localspin * nextspin * U));
} // interactions
} // basis
#pragma omp critical
total.insert(total.end(), tripletList.begin(), tripletList.end());
}
// Create the sparse Hamiltonian matrix and populate it with the triplets
observable H(num_basis, num_basis);
H.setFromTriplets(total.begin(), total.end());
return H;
}
/**
Constructs the density operator on a given site
*/
observable build_density_operator(basis basis, inversebasis inversebasis,
int site) {
int num_basis = (int)basis.size();
std::vector<T> total;
// For each basis state..
#pragma omp parallel
{
std::vector<T> tripletList;
#pragma omp for nowait
for (int i = 0; i < num_basis; ++i) {
// for (std::pair<std::string, int> element : basis) {
std::string this_state = inversebasis[i];
complex localspin = (int)(this_state[site] - '0');
tripletList.push_back(T(i, i, localspin));
}
#pragma omp critical
total.insert(total.end(), tripletList.begin(), tripletList.end());
}
// Create a new sparse matrix and populate it with the triplets
observable N(num_basis, num_basis);
N.setFromTriplets(total.begin(), total.end());
return N;
}
/**
Constructs the initial wavefunction. It is set to the
statestring specified. If that statestring does not exist,
the program will crash.
*/
state build_initial_state(basis basis, std::string statestring) {
int num_basis = (int)basis.size();
state Psi(num_basis);
Psi = Eigen::VectorXcd::Zero(num_basis);
basis::const_iterator got = basis.find(statestring);
int that_index = (int)got->second;
std::cout << " [*] State found at index: " << that_index << std::endl;
Psi.coeffRef(that_index) += complex(1, 0);
return Psi;
}
/**
* Creates an orthonormal subspace the size of Q's columns.
*/
void build_orthogonal_basis_Arnoldi(observable H, state Psi,
Eigen::MatrixXcd &Q, Eigen::MatrixXcd &h) {
Q.col(0) = Psi.normalized();
for (int i = 1; i < Q.cols(); ++i) {
// Obtain the next q vector
state q = H * Q.col(i - 1);
// Orthonormalize it w.r.t. all the others
for (int k = 0; k < i; ++k) {
h(k, i - 1) = Q.col(k).dot(q);
q = q - h(k, i - 1) * Q.col(k);
}
// Normalize it and store that in our upper Hessenberg matrix h
h(i, i - 1) = q.norm();
// Check if we found a subspace
if (h(i, i - 1).real() <= 1e-14) {
std::cerr << "Invariant subspace found before finishing Arnoldi"
<< std::endl;
}
// And set it as the next column of Q
Q.col(i) = q / h(i, i - 1);
}
// Complete h by setting the last column manually
int i = Q.cols() - 1;
h(i, i) = Q.col(i).dot(H * Q.col(i));
h(i - 1, i) = h(i, i - 1);
}
/**
* Creates an orthonormal subspace the size of Q's columns.
*/
void build_orthogonal_basis_Lanczos(observable H, state Psi,
Eigen::MatrixXcd &Q, Eigen::MatrixXcd &h) {
Eigen::VectorXcd alpha(Q.cols());
Eigen::VectorXcd beta(Q.cols());
beta(0) = 0;
Q.col(0) = Psi.normalized(); // u1; U1 = [u1]
for (int j = 1; j < Q.cols(); ++j) {
// Obtain the next q vector
state q = H * Q.col(j - 1);
alpha(j - 1) = Q.col(j - 1).dot(q);
q = q - alpha(j - 1) * Q.col(j - 1);
if (j > 1) {
q = q - beta(j - 2) * Q.col(j - 2);
}
// Re-orthogonalize
auto delta = Q.col(j - 1).dot(q);
q = q - Q.col(j - 1) * delta;
alpha(j - 1) = alpha(j - 1) + delta;
// Compute the norm
beta(j - 1) = q.norm();
// Check if we found a subspace
if (beta(j - 1).real() <= 1e-14) {
std::cerr << "Invariant subspace found before finishing Arnoldi"
<< std::endl;
}
// And set it as the next column of Q
Q.col(j) = q / beta(j - 1);
}
// Set the projected Hamiltonian
h.diagonal() = alpha.tail(Q.cols());
h.diagonal(-1) = beta.head(Q.cols() - 1);
h.diagonal(+1) = beta.head(Q.cols() - 1);
// The last diagonal has to be set separately still
int i = Q.cols() - 1;
h(i, i) = Q.col(i).dot(H * Q.col(i));
}
/**
Propagates Psi for one dt under H using a Krylov subspace
of dimension m.
*/
state krylov_propagate(observable H, state Psi, double dt, int m) {
// Krylov subspace orthogonal basis
Eigen::MatrixXcd Q(Psi.size(), m);
Q = Eigen::MatrixXcd::Zero(Psi.size(), m);
// H projected into the Krylov subspace. It is an upper Hessenberg
// matrix if doing Arnoldi, but we never need the extra bottom row
// if we're not implementing an adaptive routine that automatically
// sets m.
Eigen::MatrixXcd h(m, m);
h = Eigen::MatrixXcd::Zero(m, m);
/**
* Orthogonalize using Gram-Schmidt (sometimes called Arnoldi, since the
* to-be-orthogonalized vectors are computed on-the-fly instead of being
* stored as a matrix).
* Includes a re-orthogonalization step
* (https://link.springer.com/article/10.1007/s00211-005-0615-4)
*
* TODO: A better thing to do here might be to use householder to
* compute a QR decomposition and use that! After this first
* round of GramSchmidt-ing, the resulting vectors are not
* very orthogonal for large m (about 10).
*
* TODO: We are building a Krylov subspace of the Hamiltonian, i.e.
* of a Hermitian matrix. So we should be using Lanczos instead of
* Arnoldi!
*/
/*
// Use Arnoldi to compute Q and h
build_orthogonal_basis_Arnoldi(H, Psi, Q, h);
std::cout << "Arnoldi: " << std::endl << h << std::endl;
std::cout << "QdagQ:" << std::endl << Q.adjoint() * Q << std::endl;
*/
// Use the Lanczos algorithm with re-orthogonalization to
// compute an orthonormal Krylov subspace basis, and use it
// to obtain the projection matrix Q and the projected Hamiltonian h.
build_orthogonal_basis_Lanczos(H, Psi, Q, h);
// Convert to represent -iHdt
h *= complex(0, -dt);
// And return the first column of the back-transformed H
return (Q * h.exp()).col(0);
}
#endif //__KRYLOV_HPP__