DOCI-Exact  1.0
DM2.cpp
Go to the documentation of this file.
1 #include <algorithm>
2 #include <functional>
3 #include <hdf5.h>
4 #include <assert.h>
5 #include <omp.h>
6 #include <chrono>
7 
8 #include "DM2.h"
9 #include "DOCIHamtilonian.h"
10 #include "lapack.h"
11 
12 using namespace doci;
13 
14 std::unique_ptr<helpers::matrix> doci::DM2::sp2tp = nullptr;
15 std::unique_ptr<helpers::matrix> doci::DM2::tp2sp = nullptr;
16 
17 // this helps to check the return codes of HDF5 calls
18 #define HDF5_STATUS_CHECK(status) if(status < 0) std::cerr << __FILE__ << ":" << __LINE__ << ": Problem with writing to file. Status code=" << status << std::endl;
19 
26 DM2::DM2(unsigned int n_sp, unsigned int n)
27 {
28  if(!sp2tp || !tp2sp)
29  fill_lists(n_sp);
30 
31  block.reset(new helpers::matrix(n_sp, n_sp));
32  diag.resize((n_sp*(n_sp-1))/2);
33 
34  this->N = n;
35 }
36 
37 
43 DM2::DM2(const Molecule &mol)
44 {
45  auto n_sp = mol.get_n_sp();
46 
47  if(!sp2tp || !tp2sp)
48  fill_lists(n_sp);
49 
50  block.reset(new helpers::matrix(n_sp, n_sp));
51  diag.resize((n_sp*(n_sp-1))/2);
52  N = mol.get_n_electrons();
53 }
54 
55 DM2::DM2(const DM2 &orig)
56 {
57  assert(sp2tp && tp2sp);
58 
59  block.reset(new helpers::matrix(*orig.block));
60  diag = orig.diag;
61  N = orig.N;
62 }
63 
64 DM2::DM2(DM2 &&orig)
65 {
66  block = std::move(orig.block);
67  diag = std::move(orig.diag);
68  N = orig.N;
69 }
70 
71 DM2& DM2::operator=(const DM2 &orig)
72 {
73  block.reset(new helpers::matrix(*orig.block));
74  diag = orig.diag;
75  N = orig.N;
76 
77  return *this;
78 }
79 
81 {
82  block = std::move(orig.block);
83  diag = std::move(orig.diag);
84  N = orig.N;
85 
86  return *this;
87 }
88 
94 DM2& DM2::operator=(double val)
95 {
96  (*block) = val;
97  std::fill(diag.begin(), diag.end(), val);
98 
99  return *this;
100 }
101 
111 double DM2::operator()(int a, int b, int c, int d) const
112 {
113  if(a==b || c==d)
114  return 0;
115 
116  int sign = 1;
117 
118  if(a > b)
119  sign *= -1;
120 
121  if(c > d)
122  sign *= -1;
123 
124  int i = (*sp2tp)(a,b);
125  int j = (*sp2tp)(c,d);
126 
127  if(i<block->getn() && j<block->getn())
128  return sign * (*block)(i,j);
129  else if(i==j)
130  return sign * diag[(i-block->getn())%diag.size()];
131  else
132  return 0;
133 }
134 
141 {
142  (*block) += (*a.block);
143 // std::transform(diag.begin(), diag.end(), a.diag.begin(), diag.begin(), std::plus<double>());
144 
145  int dim = diag.size();
146  int inc = 1;
147  double alpha = 1.0;
148  double *tmp = const_cast<double *>(a.diag.data());
149 
150  daxpy_(&dim,&alpha,tmp,&inc,diag.data(),&inc);
151 
152  return *this;
153 }
154 
159 void DM2::fill_lists(unsigned int n_sp)
160 {
161  auto L = n_sp;
162  auto M = 2*n_sp;
163  auto n_tp = M*(M-1)/2;
164 
165  sp2tp.reset(new helpers::matrix(M,M));
166  (*sp2tp) = -1; // if you use something you shouldn't, this will case havoc
167 
168  tp2sp.reset(new helpers::matrix(n_tp,2));
169  (*tp2sp) = -1; // if you use something you shouldn't, this will case havoc
170 
171  auto tel = 0;
172 
173  // a \bar a
174  for(int a=0;a<L;a++)
175  (*sp2tp)(a,a+L) = (*sp2tp)(a+L,a) = tel++;
176 
177  // a b
178  for(int a=0;a<L;a++)
179  for(int b=a+1;b<L;b++)
180  (*sp2tp)(a,b) = (*sp2tp)(b,a) = tel++;
181 
182  // \bar a \bar b
183  for(int a=L;a<M;a++)
184  for(int b=a+1;b<M;b++)
185  (*sp2tp)(a,b) = (*sp2tp)(b,a) = tel++;
186 
187  // a \bar b ; a \bar b
188  for(int a=0;a<L;a++)
189  for(int b=L+a+1;b<M;b++)
190  if(a%L!=b%L)
191  (*sp2tp)(a,b) = (*sp2tp)(b,a) = tel++;
192 
193  // \bar a b ; \bar a b
194  for(int a=L;a<M;a++)
195  for(int b=a%L+1;b<L;b++)
196  if(a%L!=b%L)
197  (*sp2tp)(a,b) = (*sp2tp)(b,a) = tel++;
198 
199  assert(tel == n_tp);
200 
201  for(int a=0;a<M;a++)
202  for(int b=a+1;b<M;b++)
203  {
204  (*tp2sp)((*sp2tp)(a,b),0) = a;
205  (*tp2sp)((*sp2tp)(a,b),1) = b;
206  }
207 }
208 
213 void DM2::WriteToFile(std::string filename) const
214 {
215  hid_t file_id, group_id, dataset_id, attribute_id, dataspace_id;
216  herr_t status;
217 
218  file_id = H5Fcreate(filename.c_str(), H5F_ACC_TRUNC, H5P_DEFAULT, H5P_DEFAULT);
219  HDF5_STATUS_CHECK(file_id);
220 
221  group_id = H5Gcreate(file_id, "/RDM", H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
222  HDF5_STATUS_CHECK(group_id);
223 
224  hsize_t dimblock = block->getn() * block->getn();
225 
226  dataspace_id = H5Screate_simple(1, &dimblock, NULL);
227 
228  dataset_id = H5Dcreate(group_id, "Block", H5T_IEEE_F64LE, dataspace_id, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
229  HDF5_STATUS_CHECK(dataspace_id);
230 
231  double *data = const_cast<helpers::matrix &>(*block).getpointer();
232 
233  status = H5Dwrite(dataset_id, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, H5P_DEFAULT, data);
234  HDF5_STATUS_CHECK(status);
235 
236  status = H5Sclose(dataspace_id);
237  HDF5_STATUS_CHECK(status);
238 
239  status = H5Dclose(dataset_id);
240  HDF5_STATUS_CHECK(status);
241 
242  dimblock = diag.size();
243 
244  dataspace_id = H5Screate_simple(1, &dimblock, NULL);
245 
246  dataset_id = H5Dcreate(group_id, "Vector", H5T_IEEE_F64LE, dataspace_id, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
247  HDF5_STATUS_CHECK(dataspace_id);
248 
249  status = H5Dwrite(dataset_id, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, H5P_DEFAULT, diag.data());
250  HDF5_STATUS_CHECK(status);
251 
252  status = H5Sclose(dataspace_id);
253  HDF5_STATUS_CHECK(status);
254 
255  status = H5Dclose(dataset_id);
256  HDF5_STATUS_CHECK(status);
257 
258 
259  dataspace_id = H5Screate(H5S_SCALAR);
260  unsigned int L = block->getn();
261 
262  attribute_id = H5Acreate (group_id, "L", H5T_STD_U32LE, dataspace_id, H5P_DEFAULT, H5P_DEFAULT);
263  status = H5Awrite (attribute_id, H5T_NATIVE_UINT, &L );
264  HDF5_STATUS_CHECK(status);
265 
266  status = H5Aclose(attribute_id);
267  HDF5_STATUS_CHECK(status);
268 
269  attribute_id = H5Acreate (group_id, "N", H5T_STD_U32LE, dataspace_id, H5P_DEFAULT, H5P_DEFAULT);
270  status = H5Awrite (attribute_id, H5T_NATIVE_UINT, &N );
271  HDF5_STATUS_CHECK(status);
272 
273  status = H5Aclose(attribute_id);
274  HDF5_STATUS_CHECK(status);
275 
276  status = H5Sclose(dataspace_id);
277  HDF5_STATUS_CHECK(status);
278 
279  status = H5Gclose(group_id);
280  HDF5_STATUS_CHECK(status);
281 
282  status = H5Fclose(file_id);
283  HDF5_STATUS_CHECK(status);
284 }
285 
291 DM2 DM2::ReadFromFile(std::string filename)
292 {
293  hid_t file_id, dataset_id, group_id, attribute_id;
294  herr_t status;
295 
296  file_id = H5Fopen(filename.c_str(), H5F_ACC_RDONLY, H5P_DEFAULT);
297  HDF5_STATUS_CHECK(file_id);
298 
299  group_id = H5Gopen(file_id, "/RDM", H5P_DEFAULT);
300  HDF5_STATUS_CHECK(group_id);
301 
302  attribute_id = H5Aopen(group_id, "L", H5P_DEFAULT);
303  HDF5_STATUS_CHECK(attribute_id);
304 
305  unsigned int L;
306  status = H5Aread(attribute_id, H5T_NATIVE_UINT, &L);
307  HDF5_STATUS_CHECK(status);
308 
309  status = H5Aclose(attribute_id);
310  HDF5_STATUS_CHECK(status);
311 
312  attribute_id = H5Aopen(group_id, "N", H5P_DEFAULT);
313  HDF5_STATUS_CHECK(attribute_id);
314 
315  unsigned int N;
316  status = H5Aread(attribute_id, H5T_NATIVE_UINT, &N);
317  HDF5_STATUS_CHECK(status);
318 
319  status = H5Aclose(attribute_id);
320  HDF5_STATUS_CHECK(status);
321 
322  DM2 dm2(L,N);
323 
324  dataset_id = H5Dopen(file_id, "/RDM/Block", H5P_DEFAULT);
325  HDF5_STATUS_CHECK(dataset_id);
326 
327  status = H5Dread(dataset_id, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, H5P_DEFAULT, dm2.block->getpointer());
328  HDF5_STATUS_CHECK(status);
329 
330  status = H5Dclose(dataset_id);
331  HDF5_STATUS_CHECK(status);
332 
333  dataset_id = H5Dopen(file_id, "/RDM/Vector", H5P_DEFAULT);
334  HDF5_STATUS_CHECK(dataset_id);
335 
336  status = H5Dread(dataset_id, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, H5P_DEFAULT, dm2.diag.data());
337  HDF5_STATUS_CHECK(status);
338 
339  status = H5Dclose(dataset_id);
340  HDF5_STATUS_CHECK(status);
341 
342  status = H5Gclose(group_id);
343  HDF5_STATUS_CHECK(status);
344 
345  status = H5Fclose(file_id);
346  HDF5_STATUS_CHECK(status);
347 
348  return dm2;
349 }
350 
354 unsigned int DM2::get_n_electrons() const
355 {
356  return N;
357 }
358 
362 unsigned int DM2::get_n_sp() const
363 {
364  return block->getn();
365 }
366 
373 void DM2::Build(Permutation &perm, std::vector<double> &eigv)
374 {
375  auto num_t = omp_get_max_threads();
376  unsigned long long num_elems = (eigv.size()*1ul*(eigv.size()+1ul))/2;
377  unsigned long long size_part = num_elems/num_t + 1;
378 
379  // every thread should process the lines between i and i+1
380  // with i the thread number
381  std::vector<unsigned int> workload(num_t+1);
382  workload.front() = 0;
383  workload.back() = eigv.size();
384 
385  for(int i=1;i<num_t;i++)
386  {
387  auto num_lines = workload[i-1];
388  unsigned long long num_elems = 0;
389 
390  while(num_elems < size_part)
391  num_elems += eigv.size() - num_lines++;
392 
393  if(num_lines > eigv.size())
394  num_lines = eigv.size();
395 
396  workload[i] = num_lines;
397  }
398 
399  std::vector< std::unique_ptr<DM2> > dm2_parts(num_t);
400 
401  std::cout << "Running with " << num_t << " threads." << std::endl;
402 
403  perm.reset();
404 
405 #pragma omp parallel
406  {
407  auto start = std::chrono::high_resolution_clock::now();
408  auto me = omp_get_thread_num();
409 
410  dm2_parts[me].reset(new DM2(block->getn(),N));
411  (*dm2_parts[me]) = 0;
412 
413  Permutation my_perm(perm);
414 
415  for(auto idx_begin=0;idx_begin<workload[me];++idx_begin)
416  my_perm.next();
417 
418  auto vec_copy = eigv;
419 
420  build_iter(my_perm, vec_copy, workload[me], workload[me+1], (*dm2_parts[me]));
421 
422  auto end = std::chrono::high_resolution_clock::now();
423 
424 #pragma omp critical
425  std::cout << me << "\t" << std::chrono::duration_cast<std::chrono::duration<double,std::ratio<1>>>(end-start).count() << " s" << std::endl;
426  }
427 
428  // add everything
429  (*this) = 0;
430  for(auto &cur_dm2: dm2_parts)
431  (*this) += (*cur_dm2);
432 }
433 
434 void DM2::build_iter(Permutation& perm, std::vector<double> &eigv, unsigned int i_start, unsigned int i_end, DM2 &cur_2dm)
435 {
436  auto& perm_bra = perm;
437 
438  for(unsigned int i=i_start;i<i_end;++i)
439  {
440  const auto bra = perm_bra.get();
441 
442  Permutation perm_ket(perm_bra);
443 
444  auto cur = bra;
445 
446  // find occupied orbitals
447  while(cur)
448  {
449  // select rightmost up state in the ket
450  auto ksp = cur & (~cur + 1);
451  // set it to zero
452  cur ^= ksp;
453 
454  // number of the orbital
455  auto s = DOCIHamiltonian::CountBits(ksp-1);
456 
457  // in this case: s == (*sp2tp)(s,s+L)
458  (*cur_2dm.block)(s,s) += eigv[i] * eigv[i];
459 
460  auto cur2 = cur;
461 
462  while(cur2)
463  {
464  // select rightmost up state in the ket
465  auto ksp2 = cur2 & (~cur2 + 1);
466  // set it to zero
467  cur2 ^= ksp2;
468 
469  // number of the orbital
470  auto r = DOCIHamiltonian::CountBits(ksp2-1);
471 
472  unsigned int idx = (*sp2tp)(r,s);
473  // find correct relative index
474  idx -= block->getn();
475  idx %= diag.size();
476 
477  cur_2dm.diag[idx] += eigv[i] * eigv[i];
478  }
479  }
480 
481 
482  for(unsigned int j=i+1;j<eigv.size();++j)
483  {
484  const auto ket = perm_ket.next();
485 
486  const auto diff = bra ^ ket;
487 
488  // this means 4 orbitals are different
489  if(DOCIHamiltonian::CountBits(diff) == 2)
490  {
491  auto diff_c = diff;
492 
493  // select rightmost up state in the ket
494  auto ksp1 = diff_c & (~diff_c + 1);
495  // set it to zero
496  diff_c ^= ksp1;
497 
498  auto ksp2 = diff_c & (~diff_c + 1);
499 
500  // number of the orbital
501  auto r = DOCIHamiltonian::CountBits(ksp1-1);
502  auto s = DOCIHamiltonian::CountBits(ksp2-1);
503 
504  // in this case:
505  // r == (*sp2tp)(r,r+L)
506  // s == (*sp2tp)(s,s+L)
507  (*cur_2dm.block)(r,s) += eigv[i] * eigv[j];
508  (*cur_2dm.block)(s,r) += eigv[i] * eigv[j];
509  }
510  }
511 
512  perm_bra.next();
513  }
514 }
515 
516 std::ostream &operator<<(std::ostream &output,doci::DM2 &dm2)
517 {
518  auto L = dm2.get_n_sp();
519  output << "Block: " << std::endl;
520  for(int i=0;i<L;i++)
521  for(int j=i;j<L;j++)
522  output << i << "\t" << j << "\t|\t" << (*dm2.tp2sp)(i,0) << " " << (*dm2.tp2sp)(i,1) << " ; " << (*dm2.tp2sp)(j,0) << " " << (*dm2.tp2sp)(j,1) << "\t\t" << (*dm2.block)(i,j) << std::endl;
523 
524  output << std::endl;
525 
526  output << "Vector (4x): " << std::endl;
527  for(int i=0;i<dm2.diag.size();i++)
528  output << i << "\t|\t" << (*dm2.tp2sp)(L+i,0) << " " << (*dm2.tp2sp)(L+i,1) << "\t\t" << dm2.diag[i] << std::endl;
529 
530  return output;
531 }
532 
538 {
539  auto L = block->getn();
540 
541  // make our life easier
542  auto calc_elem = [this,&mol,L] (int i, int j) {
543  int a = (*tp2sp)(i,0);
544  int b = (*tp2sp)(i,1);
545  int c = (*tp2sp)(j,0);
546  int d = (*tp2sp)(j,1);
547 
548  int a_ = a % L;
549  int b_ = b % L;
550  int c_ = c % L;
551  int d_ = d % L;
552 
553  double result = 0;
554 
555  // sp terms
556  if(i==j)
557  result += (mol.getT(a_,a_) + mol.getT(b_,b_))/(N - 1.0);
558 
559  // tp terms:
560 
561  // a \bar a ; b \bar b
562  if(b==(a+L) && d==(c+L))
563  result += mol.getV(a_,b_,c_,d_);
564 
565  // a b ; a b
566  // \bar a \bar b ; \bar a \bar b
567  if(i==j && a/L == b/L && a!=b)
568  result += mol.getV(a_,b_,c_,d_) - mol.getV(a_,b_,d_,c_);
569 
570  // a \bar b ; a \bar b
571  // \bar a b ; \bar a b
572  if(i==j && a/L != b/L && a%L!=b%L)
573  result += mol.getV(a_,b_, c_,d_);
574 
575  return result;
576  };
577 
578  for(int i=0;i<block->getn();++i)
579  for(int j=i;j<block->getn();++j)
580  (*block)(i,j) = (*block)(j,i) = calc_elem(i,j);
581 
582  for(int i=0;i<diag.size();i++)
583  // keep in mind that the degen of the vector is 4. We need prefactor of 2, so
584  // we end up with 0.5
585  diag[i] = 0.5*calc_elem(L+i,L+i) + 0.5*calc_elem(L*L+i,L*L+i);
586 }
587 
593 double DM2::Dot(const DM2 &x) const
594 {
595  double result = 0;
596 
597  int n = block->getn()*block->getn();
598  int inc = 1;
599 
600  result += ddot_(&n,block->getpointer(),&inc,x.block->getpointer(),&inc);
601 
602  n = diag.size();
603  result += 4 * ddot_(&n,diag.data(),&inc,x.diag.data(),&inc);
604 
605  return result;
606 }
607 
611 double DM2::Trace() const
612 {
613  double result = 0;
614  for(auto elem : diag)
615  result += elem;
616 
617  return block->trace() + 4 * result;
618 }
619 
630 std::pair<double,bool> DM2::find_min_angle(int k, int l, double start_angle, std::function<double(int,int)> &T, std::function<double(int,int,int,int)> &V) const
631 {
632  assert(k!=l);
633 
634  const int L = block->getn();
635 
636  double theta = start_angle;
637 
638  const DM2 &rdm = *this;
639 
640  double cos2 = 2.0/(N-1.0)*(T(k,k)*rdm(k,k+L,k,k+L)+T(l,l)*rdm(l,l+L,l,l+L));
641 
642  double sin2 = 2.0/(N-1.0)*(T(l,l)*rdm(k,k+L,k,k+L)+T(k,k)*rdm(l,l+L,l,l+L));
643 
644  // 2sincos actually
645  double sincos = 2.0/(N-1.0)*T(k,l)*(rdm(l,l+L,l,l+L)-rdm(k,k+L,k,k+L));
646 
647  for(int a=0;a<L;a++)
648  {
649  if(a==k || a==l)
650  continue;
651 
652  cos2 += 2*V(k,k,a,a)*rdm(k,k+L,a,a+L)+2*V(l,l,a,a)*rdm(l,l+L,a,a+L)+2*(2*V(k,a,k,a)-V(k,a,a,k)+2.0/(N-1.0)*T(k,k))*rdm(k,a,k,a)+2*(2*V(l,a,l,a)-V(l,a,a,l)+2.0/(N-1.0)*T(l,l))*rdm(l,a,l,a);
653 
654  sin2 += 2*V(l,l,a,a)*rdm(k,k+L,a,a+L)+2*V(k,k,a,a)*rdm(l,l+L,a,a+L)+2*(2*V(k,a,k,a)-V(k,a,a,k)+2.0/(N-1.0)*T(k,k))*rdm(l,a,l,a)+2*(2*V(l,a,l,a)-V(l,a,a,l)+2.0/(N-1.0)*T(l,l))*rdm(k,a,k,a);
655 
656  sincos += 2*V(k,l,a,a)*(rdm(l,l+L,a,a+L)-rdm(k,k+L,a,a+L))+2*(2*V(k,a,l,a)-V(k,a,a,l)+2.0/(N-1.0)*T(k,l))*(rdm(l,a,l,a)-rdm(k,a,k,a));
657  }
658 
659  const double cos4 = V(k,k,k,k)*rdm(k,k+L,k,k+L)+V(l,l,l,l)*rdm(l,l+L,l,l+L)+2*V(k,k,l,l)*rdm(k,k+L,l,l+L)+2*(2*V(k,l,k,l)-V(k,k,l,l))*rdm(k,l,k,l);
660 
661  const double sin4 = V(k,k,k,k)*rdm(l,l+L,l,l+L)+V(l,l,l,l)*rdm(k,k+L,k,k+L)+2*V(k,k,l,l)*rdm(k,k+L,l,l+L)+2*(2*V(k,l,k,l)-V(k,k,l,l))*rdm(k,l,k,l);
662 
663  // 2 x
664  const double cos2sin2 = (2*V(k,k,l,l)+V(k,l,k,l))*(rdm(k,k+L,k,k+L)+rdm(l,l+L,l,l+L))+((V(k,k,k,k)+V(l,l,l,l)-2*(V(k,l,k,l)+V(k,k,l,l))))*rdm(k,k+L,l,l+L)+(V(k,k,k,k)+V(l,l,l,l)-6*V(k,k,l,l)+2*V(k,l,k,l))*rdm(k,l,k,l);
665 
666  // 4 x
667  const double sin3cos = V(k,l,k,k)*rdm(l,l+L,l,l+L)-V(k,l,l,l)*rdm(k,k+L,k,k+L)-(V(k,l,k,k)-V(k,l,l,l))*(rdm(k,k+L,l,l+L)+rdm(k,l,k,l));
668 
669  // 4 x
670  const double cos3sin = V(k,l,l,l)*rdm(l,l+L,l,l+L)-V(k,l,k,k)*rdm(k,k+L,k,k+L)+(V(k,l,k,k)-V(k,l,l,l))*(rdm(k,k+L,l,l+L)+rdm(k,l,k,l));
671 
672  // A*cos(t)^4+B*sin(t)^4+C*cos(t)^2+D*sin(t)^2+2*E*cos(t)*sin(t)+2*F*cos(t)^2*sin(t)^2+4*G*sin(t)*cos(t)^3+4*H*sin(t)^3*cos(t)
673 
674  // (16*G-16*H)*cos(t)^4+(-4*A-4*B+8*F)*sin(t)*cos(t)^3+(4*E-12*G+20*H)*cos(t)^2+(4*B-2*C+2*D-4*F)*sin(t)*cos(t)-2*E-4*H
675  auto gradient = [&] (double theta) -> double {
676  double cos = std::cos(theta);
677  double sin = std::sin(theta);
678 
679  double res = 16*(cos3sin-sin3cos)*cos*cos*cos*cos - 4*(cos4+sin4-2*cos2sin2)*sin*cos*cos*cos + 4*(sincos-3*cos3sin+5*sin3cos)*cos*cos + 2*(2*sin4-cos2+sin2-2*cos2sin2)*sin*cos-2*sincos-4*sin3cos;
680 
681  return res;
682  };
683 
684  // (-16*A-16*B+32*F)*cos(t)^4+(-64*G+64*H)*sin(t)*cos(t)^3+(12*A+20*B-4*C+4*D-32*F)*cos(t)^2+(-8*E+24*G-40*H)*sin(t)*cos(t)-4*B+2*C-2*D+4*F
685  auto hessian = [&] (double theta) -> double {
686  double cos = std::cos(theta);
687  double sin = std::sin(theta);
688 
689  double res = -16*(cos4+sin4-2*cos2sin2)*cos*cos*cos*cos+64*(sin3cos-cos3sin)*sin*cos*cos*cos+4*(3*cos4+5*sin4-cos2+sin2-8*cos2sin2)*cos*cos-8*(sincos-3*cos3sin+5*sin3cos)*sin*cos+2*(cos2-2*sin4-sin2+2*cos2sin2);
690 
691  return res;
692  };
693 
694 // int Na = 5000;
695 // for(int i=0;i<Na;i++)
696 // {
697 // double t = 2.0 * M_PI / (1.0*Na) * i;
698 // std::cout << t << "\t" << calc_rotate(k,l,t,T,V) << "\t" << gradient(t) << "\t" << hessian(t) << std::endl;
699 // }
700 
701  const int max_iters = 20;
702  const double convergence = 1e-12;
703 
704  double change = gradient(theta)*theta+hessian(theta)*theta*theta/2.0;
705 
706  // if it goes uphill, try flipping sign and try again
707  if(change>0)
708  theta *= -1;
709 
710  int iter;
711  for(iter=0;iter<max_iters;iter++)
712  {
713  double dx = gradient(theta)/hessian(theta);
714 
715  theta -= dx;
716 
717  if(fabs(dx) < convergence)
718  break;
719  }
720 
721  if(iter>=max_iters)
722  std::cout << "Reached max iters in find_min_angle: " << k << " " << l << std::endl;
723 
724  if(hessian(theta)<0)
725  std::cout << "Found max!" << std::endl;
726 
727  return std::make_pair(theta, hessian(theta)>0);
728 }
729 
740 double DM2::calc_rotate(int k, int l, double theta, std::function<double(int,int)> &T, std::function<double(int,int,int,int)> &V) const
741 {
742  assert(k!=l);
743 
744  const int L = block->getn();
745 
746  const DM2 &rdm = *this;
747 
748  double energy = 4/(N-1.0)*(T(k,k)+T(l,l)) * rdm(k,l,k,l);
749 
750  double cos2 = 2.0/(N-1.0)*(T(k,k)*rdm(k,k+L,k,k+L)+T(l,l)*rdm(l,l+L,l,l+L));
751 
752  double sin2 = 2.0/(N-1.0)*(T(l,l)*rdm(k,k+L,k,k+L)+T(k,k)*rdm(l,l+L,l,l+L));
753 
754  // 2sincos actually
755  double sincos = 2.0/(N-1.0)*T(k,l)*(rdm(l,l+L,l,l+L)-rdm(k,k+L,k,k+L));
756 
757  for(int a=0;a<L;a++)
758  {
759  if(a==k || a==l)
760  continue;
761 
762  energy += 2.0/(N-1.0) * T(a,a) * (rdm(a,a+L,a,a+L)+2*rdm(a,k,a,k)+2*rdm(a,l,a,l));
763 
764  for(int b=0;b<L;b++)
765  {
766  if(b==k || b==l)
767  continue;
768 
769  energy += 2.0/(N-1.0) * (T(a,a)+T(b,b)) * rdm(a,b,a,b);
770 
771  energy += V(a,a,b,b) * rdm(a,a+L,b,b+L);
772 
773  energy += (2*V(a,b,a,b)-V(a,b,b,a)) * rdm(a,b,a,b);
774  }
775 
776  cos2 += 2*V(k,k,a,a)*rdm(k,k+L,a,a+L)+2*V(l,l,a,a)*rdm(l,l+L,a,a+L)+2*(2*V(k,a,k,a)-V(k,a,a,k)+2.0/(N-1.0)*T(k,k))*rdm(k,a,k,a)+2*(2*V(l,a,l,a)-V(l,a,a,l)+2.0/(N-1.0)*T(l,l))*rdm(l,a,l,a);
777 
778  sin2 += 2*V(l,l,a,a)*rdm(k,k+L,a,a+L)+2*V(k,k,a,a)*rdm(l,l+L,a,a+L)+2*(2*V(k,a,k,a)-V(k,a,a,k)+2.0/(N-1.0)*T(k,k))*rdm(l,a,l,a)+2*(2*V(l,a,l,a)-V(l,a,a,l)+2.0/(N-1.0)*T(l,l))*rdm(k,a,k,a);
779 
780  sincos += 2*V(k,l,a,a)*(rdm(l,l+L,a,a+L)-rdm(k,k+L,a,a+L))+2*(2*V(k,a,l,a)-V(k,a,a,l)+2.0/(N-1.0)*T(k,l))*(rdm(l,a,l,a)-rdm(k,a,k,a));
781  }
782 
783  const double cos4 = V(k,k,k,k)*rdm(k,k+L,k,k+L)+V(l,l,l,l)*rdm(l,l+L,l,l+L)+2*V(k,k,l,l)*rdm(k,k+L,l,l+L)+2*(2*V(k,l,k,l)-V(k,k,l,l))*rdm(k,l,k,l);
784 
785  const double sin4 = V(k,k,k,k)*rdm(l,l+L,l,l+L)+V(l,l,l,l)*rdm(k,k+L,k,k+L)+2*V(k,k,l,l)*rdm(k,k+L,l,l+L)+2*(2*V(k,l,k,l)-V(k,k,l,l))*rdm(k,l,k,l);
786 
787  // 2 x
788  const double cos2sin2 = (2*V(k,k,l,l)+V(k,l,k,l))*(rdm(k,k+L,k,k+L)+rdm(l,l+L,l,l+L))+((V(k,k,k,k)+V(l,l,l,l)-2*(V(k,l,k,l)+V(k,k,l,l))))*rdm(k,k+L,l,l+L)+(V(k,k,k,k)+V(l,l,l,l)-6*V(k,k,l,l)+2*V(k,l,k,l))*rdm(k,l,k,l);
789 
790  // 4 x
791  const double sin3cos = V(k,l,k,k)*rdm(l,l+L,l,l+L)-V(k,l,l,l)*rdm(k,k+L,k,k+L)-(V(k,l,k,k)-V(k,l,l,l))*(rdm(k,k+L,l,l+L)+rdm(k,l,k,l));
792 
793  // 4 x
794  const double cos3sin = V(k,l,l,l)*rdm(l,l+L,l,l+L)-V(k,l,k,k)*rdm(k,k+L,k,k+L)+(V(k,l,k,k)-V(k,l,l,l))*(rdm(k,k+L,l,l+L)+rdm(k,l,k,l));
795 
796  const double cos = std::cos(theta);
797  const double sin = std::sin(theta);
798 
799  energy += cos*cos*cos*cos*cos4;
800 
801  energy += sin*sin*sin*sin*sin4;
802 
803  energy += cos*cos*cos2;
804 
805  energy += sin*sin*sin2;
806 
807  energy += 2*sin*cos*sincos;
808 
809  energy += 2*cos*cos*sin*sin*cos2sin2;
810 
811  energy += 4*cos*sin*sin*sin*sin3cos;
812 
813  energy += 4*cos*cos*cos*sin*cos3sin;
814 
815  return energy;
816 }
817 
818 /* vim: set ts=3 sw=3 expandtab :*/
virtual mybitset next()
Definition: Permutation.cpp:30
std::pair< double, bool > find_min_angle(int k, int l, double start_angle, std::function< double(int, int)> &T, std::function< double(int, int, int, int)> &V) const
Definition: DM2.cpp:630
static DM2 ReadFromFile(const std::string)
Definition: DM2.cpp:291
double Trace() const
Definition: DM2.cpp:611
double calc_rotate(int k, int l, double theta, std::function< double(int, int)> &T, std::function< double(int, int, int, int)> &V) const
Definition: DM2.cpp:740
double ddot_(int *n, double *x, int *incx, double *y, int *incy)
void WriteToFile(const std::string) const
Definition: DM2.cpp:213
double operator()(int, int, int, int) const
Definition: DM2.cpp:111
Definition: DM2.h:24
static unsigned int CountBits(mybitset)
unsigned int get_n_sp() const
Definition: DM2.cpp:362
std::ostream & operator<<(std::ostream &output, doci::DM2 &dm2)
Definition: DM2.cpp:516
virtual mybitset get() const
Definition: Permutation.cpp:57
void Build(Permutation &, std::vector< double > &)
Definition: DM2.cpp:373
void BuildHamiltonian(const Molecule &)
Definition: DM2.cpp:537
DM2(unsigned int, unsigned int)
Definition: DM2.cpp:26
virtual unsigned int get_n_electrons() const
Definition: Molecule.cpp:32
virtual double getT(int, int) const =0
void daxpy_(int *n, double *alpha, double *x, int *incx, double *y, int *incy)
DM2 & operator=(const DM2 &)
Definition: DM2.cpp:71
#define HDF5_STATUS_CHECK(status)
Definition: DM2.cpp:18
virtual double getV(int, int, int, int) const =0
unsigned int get_n_electrons() const
Definition: DM2.cpp:354
virtual void reset()
Definition: Permutation.cpp:66
double Dot(const DM2 &) const
Definition: DM2.cpp:593
double * getpointer() const
Definition: helpers.cpp:165
DM2 & operator+=(const DM2 &)
Definition: DM2.cpp:140
virtual unsigned int get_n_sp() const
Definition: Molecule.cpp:24