HubbardGPU
Hubbard diagonalisation on the GPU (and CPU)
 All Classes Files Functions Variables Typedefs Friends Macros
ham-mom.cpp
Go to the documentation of this file.
1 /* Copyright (C) 2014 Ward Poelmans
2 
3 This file is part of Hubbard-GPU.
4 
5 Hubbard-GPU is free software: you can redistribute it and/or modify
6 it under the terms of the GNU General Public License as published by
7 the Free Software Foundation, either version 3 of the License, or
8 (at your option) any later version.
9 
10 Hubbard-GPU is distributed in the hope that it will be useful,
11 but WITHOUT ANY WARRANTY; without even the implied warranty of
12 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 GNU General Public License for more details.
14 
15 You should have received a copy of the GNU General Public License
16 along with Hubbard-GPU. If not, see <http://www.gnu.org/licenses/>.
17 */
18 
19 #include <iostream>
20 #include <algorithm>
21 #include <cmath>
22 #include <sstream>
23 #include <cstring>
24 #include <tuple>
25 #include <iomanip>
26 #include <hdf5.h>
27 #include "ham-mom.h"
28 
29 // macro to help check return status of HDF5 functions
30 #define HDF5_STATUS_CHECK(status) { \
31  if(status < 0) \
32  std::cerr << __FILE__ << ":" << __LINE__ << \
33  ": Problem with writing to file. Status code=" \
34  << status << std::endl; \
35 }
36 
45 MomHamiltonian::MomHamiltonian(int L, int Nu, int Nd, double J, double U): BareHamiltonian(L,Nu,Nd,J,U)
46 {
47 }
48 
53 {
54 }
55 
61 {
62  baseUp.reserve(CalcDim(L,Nu));
63  baseDown.reserve(CalcDim(L,Nd));
64 
65  for(myint i=0;i<Hb;i++)
66  {
67  if(CountBits(i) == Nd)
68  baseDown.push_back(i);
69 
70  if(CountBits(i) == Nu)
71  baseUp.push_back(i);
72  }
73 
74  std::vector< std::tuple<int,int,int> > totalmom;
75  totalmom.reserve(dim);
76 
77  // count momentum of earch state
78  for(unsigned int a=0;a<baseUp.size();a++)
79  for(unsigned int b=0;b<baseDown.size();b++)
80  {
81  auto calcK = [] (myint cur)
82  {
83  int tot = 0;
84 
85  while(cur)
86  {
87  // select rightmost up state in the ket
88  myint ksp = cur & (~cur + 1);
89  // set it to zero
90  cur ^= ksp;
91 
92  tot += BareHamiltonian::CountBits(ksp-1);
93  }
94 
95  return tot;
96  };
97 
98  int K = calcK(baseUp[a]) + calcK(baseDown[b]);
99 
100  K = K % L;
101 
102  totalmom.push_back(std::make_tuple(a,b,K));
103  }
104 
105  std::sort(totalmom.begin(), totalmom.end(),
106  [](const std::tuple<int,int,int> & a, const std::tuple<int,int,int> & b) -> bool
107  {
108  return std::get<2>(a) < std::get<2>(b);
109  });
110 
111  // a block for each momenta
112  mombase.resize(L);
113 
114  std::for_each(totalmom.begin(), totalmom.end(), [this](std::tuple<int,int,int> elem)
115  {
116  auto tmp = std::make_pair(std::get<0>(elem), std::get<1>(elem));
117  mombase[std::get<2>(elem)].push_back(tmp);
118  } );
119 
120 // for(unsigned int i=0;i<mombase.size();i++)
121 // {
122 // std::cout << "K = " << i << " (" << mombase[i].size() << ")" << std::endl;
123 // for(unsigned int j=0;j<mombase[i].size();j++)
124 // std::cout << j << "\t" << mombase[i][j].first << "\t" << mombase[i][j].second << "\t" << print_bin(baseUp[mombase[i][j].first]) << "\t" << print_bin(baseDown[mombase[i][j].second]) << std::endl;
125 // }
126 }
127 
132 {
133  if( !baseUp.size() || !baseDown.size() )
134  {
135  std::cerr << "Build base before building MomHamiltonian" << std::endl;
136  return;
137  }
138 
139  blockmat.resize(L);
140 
141  for(int B=0;B<L;B++)
142  {
143  int cur_dim = mombase[B].size();
144 
145  blockmat[B].reset(new double [cur_dim*cur_dim]);
146 
147  for(int i=0;i<cur_dim;i++)
148  {
149  int a = mombase[B][i].first;
150  int b = mombase[B][i].second;
151 
152  for(int j=i;j<cur_dim;j++)
153  {
154  int c = mombase[B][j].first;
155  int d = mombase[B][j].second;
156 
157  blockmat[B][j+cur_dim*i] = U/L*interaction(a,b,c,d);
158 
159  blockmat[B][i+cur_dim*j] = blockmat[B][j+cur_dim*i];
160 
161  }
162 
163  blockmat[B][i+cur_dim*i] += J * (hopping(baseUp[a]) + hopping(baseDown[b]));
164  }
165  }
166 }
167 
177 int MomHamiltonian::interaction(int a, int b, int c, int d) const
178 {
179  int result = 0;
180 
181  myint downket = baseDown[d];
182  myint upket, upbra;
183 
184  while(downket)
185  {
186  // select rightmost down state in the ket
187  myint k2sp = downket & (~downket + 1);
188  // set it to zero
189  downket ^= k2sp;
190 
191  int k2 = CountBits(k2sp-1);
192 
193  int signk2 = 1;
194 
195  // when uneven, extra minus sign
196  if(Nu & 0x1)
197  signk2 *= -1;
198 
199  signk2 *= CalcSign(k2,Hbc,baseDown[d]);
200 
201  upket = baseUp[c];
202 
203  while(upket)
204  {
205  // select rightmost up state in the ket
206  myint k1sp = upket & (~upket + 1);
207  // set it to zero
208  upket ^= k1sp;
209 
210  int k1 = CountBits(k1sp-1);
211 
212  int signk1 = CalcSign(k1,Hbc,baseUp[c]);
213 
214  myint K = (k1+k2) % L;
215 
216  upbra = baseUp[a];
217 
218  while(upbra)
219  {
220  // select rightmost up state in the ket
221  myint k3sp = upbra & (~upbra + 1);
222  // set it to zero
223  upbra ^= k3sp;
224 
225  // k3: spin up in bra
226  int k3 = CountBits(k3sp-1);
227 
228  if(! ((baseUp[c] ^ k1sp) == (baseUp[a] ^ k3sp)) )
229  continue;
230 
231  // k4: spin down in bra
232  int k4 = (K-k3+L)%L;
233 
234  if(! ((1<<k4 & baseDown[b]) && ((baseDown[d] ^ k2sp) == (baseDown[b] ^ 1<<k4)) ))
235  continue;
236 
237  int signk3 = CalcSign(k3,Hbc,baseUp[a]);
238 
239  int signk4 = 1;
240 
241  // when uneven, extra minus sign
242  if(Nu & 0x1)
243  signk4 *= -1;
244 
245  signk4 *= CalcSign(k4,Hbc,baseDown[b]);
246 
247  result += signk1*signk2*signk3*signk4;
248  }
249  }
250  }
251 
252  return result;
253 }
254 
263 {
264  double result = 0;
265 
266  while(ket)
267  {
268  // select rightmost one
269  myint hop = ket & (~ket + 1);
270  // set it to zero in ket
271  ket ^= hop;
272 
273  result += -2* cos(2.0*M_PI/L * CountBits(hop-1));
274  }
275 
276  return result;
277 }
278 
285 void MomHamiltonian::mvprod(double *x, double *y, double alpha) const
286 {
287  double beta = 1;
288  int incx = 1;
289  char uplo = 'U';
290 
291  int offset = 0;
292 
293  for(int B=0;B<L;B++)
294  {
295  int dim = mombase[B].size();
296 
297  dsymv_(&uplo,&dim,&beta,blockmat[B].get(),&dim,&x[offset],&incx,&alpha,&y[offset],&incx);
298 
299  offset += dim;
300  }
301 }
302 
310 std::vector<double> MomHamiltonian::ExactDiagonalizeFull(bool calc_eigenvectors)
311 {
312  std::vector<double> eigenvalues(dim);
313 
314  int offset = 0;
315 
316  for(int B=0;B<L;B++)
317  {
318  int dim = mombase[B].size();
319 
320  Diagonalize(dim, blockmat[B].get(), &eigenvalues[offset], calc_eigenvectors);
321 
322  offset += dim;
323  }
324 
325  std::sort(eigenvalues.begin(), eigenvalues.end());
326 
327  return eigenvalues;
328 }
329 
338 std::vector< std::pair<int,double> > MomHamiltonian::ExactMomDiagonalizeFull(bool calc_eigenvectors)
339 {
340  std::vector<double> eigs(dim);
341  std::vector< std::pair<int, double> > energy;
342 
343  int offset = 0;
344 
345  for(int B=0;B<L;B++)
346  {
347  int dim = mombase[B].size();
348 
349  Diagonalize(dim, blockmat[B].get(), &eigs[offset], calc_eigenvectors);
350 
351  for(int i=0;i<dim;i++)
352  energy.push_back(std::make_pair(B,eigs[offset+i]));
353 
354  offset += dim;
355  }
356 
357  std::sort(energy.begin(), energy.end(),
358  [](const std::pair<int,double> & a, const std::pair<int,double> & b) -> bool
359  {
360  return a.second < b.second;
361  });
362 
363  return energy;
364 }
365 
374 void MomHamiltonian::GenerateData(double Ubegin, double Uend, double step, std::string filename)
375 {
376  double Ucur = Ubegin;
377 
378  if( !baseUp.size() || !baseDown.size() )
379  BuildBase();
380 
381  std::vector<double> eigenvalues(dim);
382 
383  hid_t file_id, group_id, dataset_id, dataspace_id, attribute_id;
384  herr_t status;
385 
386  file_id = H5Fcreate(filename.c_str(), H5F_ACC_TRUNC, H5P_DEFAULT, H5P_DEFAULT);
387  HDF5_STATUS_CHECK(file_id);
388 
389  group_id = H5Gcreate(file_id, "run", H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
390  HDF5_STATUS_CHECK(group_id);
391 
392  dataspace_id = H5Screate(H5S_SCALAR);
393 
394  attribute_id = H5Acreate (group_id, "L", H5T_STD_I32LE, dataspace_id, H5P_DEFAULT, H5P_DEFAULT);
395  status = H5Awrite (attribute_id, H5T_NATIVE_INT, &L );
396  HDF5_STATUS_CHECK(status);
397  status = H5Aclose(attribute_id);
398  HDF5_STATUS_CHECK(status);
399 
400  attribute_id = H5Acreate (group_id, "Nu", H5T_STD_I32LE, dataspace_id, H5P_DEFAULT, H5P_DEFAULT);
401  status = H5Awrite (attribute_id, H5T_NATIVE_INT, &Nu );
402  HDF5_STATUS_CHECK(status);
403  status = H5Aclose(attribute_id);
404  HDF5_STATUS_CHECK(status);
405 
406  attribute_id = H5Acreate (group_id, "Nd", H5T_STD_I32LE, dataspace_id, H5P_DEFAULT, H5P_DEFAULT);
407  status = H5Awrite (attribute_id, H5T_NATIVE_INT, &Nd );
408  HDF5_STATUS_CHECK(status);
409  status = H5Aclose(attribute_id);
410  HDF5_STATUS_CHECK(status);
411 
412  attribute_id = H5Acreate (group_id, "J", H5T_IEEE_F64LE, dataspace_id, H5P_DEFAULT, H5P_DEFAULT);
413  status = H5Awrite (attribute_id, H5T_NATIVE_DOUBLE, &J );
414  HDF5_STATUS_CHECK(status);
415  status = H5Aclose(attribute_id);
416  HDF5_STATUS_CHECK(status);
417 
418  attribute_id = H5Acreate (group_id, "Ubegin", H5T_IEEE_F64LE, dataspace_id, H5P_DEFAULT, H5P_DEFAULT);
419  status = H5Awrite (attribute_id, H5T_NATIVE_DOUBLE, &Ubegin );
420  HDF5_STATUS_CHECK(status);
421  status = H5Aclose(attribute_id);
422  HDF5_STATUS_CHECK(status);
423 
424  attribute_id = H5Acreate (group_id, "Uend", H5T_IEEE_F64LE, dataspace_id, H5P_DEFAULT, H5P_DEFAULT);
425  status = H5Awrite (attribute_id, H5T_NATIVE_DOUBLE, &Uend );
426  HDF5_STATUS_CHECK(status);
427  status = H5Aclose(attribute_id);
428  HDF5_STATUS_CHECK(status);
429 
430  attribute_id = H5Acreate (group_id, "Ustep", H5T_IEEE_F64LE, dataspace_id, H5P_DEFAULT, H5P_DEFAULT);
431  status = H5Awrite (attribute_id, H5T_NATIVE_DOUBLE, &step );
432  HDF5_STATUS_CHECK(status);
433  status = H5Aclose(attribute_id);
434  HDF5_STATUS_CHECK(status);
435 
436  status = H5Sclose(dataspace_id);
437  HDF5_STATUS_CHECK(status);
438 
439  status = H5Gclose(group_id);
440  HDF5_STATUS_CHECK(status);
441 
442  status = H5Fclose(file_id);
443  HDF5_STATUS_CHECK(status);
444 
445  std::vector<double> diagonalelements(dim);
446 
447  std::vector< std::unique_ptr<double []> > offdiag;
448  offdiag.resize(L);
449 
450  // make sure that we don't rebuild the whole hamiltonian every time.
451  // store the hopping and interaction part seperate so we can just
452  // add them in every step
453 #pragma omp parallel for
454  for(int B=0;B<L;B++)
455  {
456  int cur_dim = mombase[B].size();
457  int offset = 0;
458 
459  for(int tmp=0;tmp<B;tmp++)
460  offset += mombase[tmp].size();
461 
462  offdiag[B].reset(new double [cur_dim*cur_dim]);
463 
464  for(int i=0;i<cur_dim;i++)
465  {
466  int a = mombase[B][i].first;
467  int b = mombase[B][i].second;
468 
469  diagonalelements[offset+i] = hopping(baseUp[a]) + hopping(baseDown[b]);
470 
471  for(int j=i;j<cur_dim;j++)
472  {
473  int c = mombase[B][j].first;
474  int d = mombase[B][j].second;
475 
476  offdiag[B][j+cur_dim*i] = 1.0/L*interaction(a,b,c,d);
477  offdiag[B][i+cur_dim*j] = offdiag[B][j+cur_dim*i];
478  }
479  }
480  }
481 
482 
483  while(Ucur <= Uend)
484  {
485  std::cout << "U = " << Ucur << std::endl;
486  setU(Ucur);
487 
488  // make hamiltonian
489 #pragma omp parallel for
490  for(int B=0;B<L;B++)
491  {
492  int cur_dim = mombase[B].size();
493  int offset = 0;
494 
495  for(int tmp=0;tmp<B;tmp++)
496  offset += mombase[tmp].size();
497 
498  std::memcpy(blockmat[B].get(),offdiag[B].get(),cur_dim*cur_dim*sizeof(double));
499 
500  int tmp = cur_dim*cur_dim;
501  int inc = 1;
502  dscal_(&tmp,&Ucur,blockmat[B].get(),&inc);
503 
504  for(int i=0;i<cur_dim;i++)
505  blockmat[B][i+cur_dim*i] += diagonalelements[offset+i];
506  }
507 
508 
509 #pragma omp parallel for
510  for(int B=0;B<L;B++)
511  {
512  int dim = mombase[B].size();
513  int offset = 0;
514 
515  for(int tmp=0;tmp<B;tmp++)
516  offset += mombase[tmp].size();
517 
518  Diagonalize(dim, blockmat[B].get(), &eigenvalues[offset], false);
519  }
520 
521  hid_t U_id;
522  std::stringstream name;
523  name << std::setprecision(5) << std::fixed << "/run/" << Ucur;
524 
525  file_id = H5Fopen(filename.c_str(), H5F_ACC_RDWR, H5P_DEFAULT);
526  HDF5_STATUS_CHECK(file_id);
527 
528  U_id = H5Gcreate(file_id, name.str().c_str(), H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
529  HDF5_STATUS_CHECK(U_id);
530 
531  for(int B=0;B<L;B++)
532  {
533  int dim = mombase[B].size();
534  int offset = 0;
535 
536  for(int tmp=0;tmp<B;tmp++)
537  offset += mombase[tmp].size();
538 
539  hsize_t dimarr = dim;
540 
541  dataspace_id = H5Screate_simple(1, &dimarr, NULL);
542 
543  std::stringstream cur_block;
544  cur_block << B;
545  dataset_id = H5Dcreate(U_id, cur_block.str().c_str(), H5T_IEEE_F64LE, dataspace_id, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
546 
547  status = H5Dwrite(dataset_id, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, H5P_DEFAULT, &eigenvalues[offset] );
548  HDF5_STATUS_CHECK(status);
549 
550  status = H5Sclose(dataspace_id);
551  HDF5_STATUS_CHECK(status);
552 
553  status = H5Dclose(dataset_id);
554  HDF5_STATUS_CHECK(status);
555  }
556 
557  status = H5Gclose(U_id);
558  HDF5_STATUS_CHECK(status);
559 
560  status = H5Fclose(file_id);
561  HDF5_STATUS_CHECK(status);
562 
563  Ucur += step;
564  }
565 }
566 
567 /* vim: set ts=8 sw=4 tw=0 expandtab :*/