CheMPS2
DMRGtechnics.cpp
1 /*
2  CheMPS2: a spin-adapted implementation of DMRG for ab initio quantum chemistry
3  Copyright (C) 2013-2016 Sebastian Wouters
4 
5  This program 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 2 of the License, or
8  (at your option) any later version.
9 
10  This program 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 along
16  with this program; if not, write to the Free Software Foundation, Inc.,
17  51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
18 */
19 
20 #include <iostream>
21 #include <stdlib.h>
22 #include <sstream>
23 #include <string>
24 #include <algorithm>
25 #include <math.h>
26 #include <assert.h>
27 #include <sys/time.h>
28 #include <unistd.h>
29 
30 #include "DMRG.h"
31 #include "Lapack.h"
32 #include "Heff.h"
33 #include "MPIchemps2.h"
34 #include "Wigner.h"
35 #include "Special.h"
36 
37 using std::cout;
38 using std::endl;
39 
40 void CheMPS2::DMRG::calc_rdms_and_correlations( const bool do_3rdm, const bool disk_3rdm ){
41 
42  #ifdef CHEMPS2_MPI_COMPILATION
43  const bool am_i_master = ( MPIchemps2::mpi_rank() == MPI_CHEMPS2_MASTER );
44  #else
45  const bool am_i_master = true;
46  #endif
47 
48  // Reset timings
49  for ( int timecnt = 0; timecnt < CHEMPS2_TIME_VECLENGTH; timecnt++ ){ timings[ timecnt ] = 0.0; }
50  num_double_write_disk = 0;
51  num_double_read_disk = 0;
52  struct timeval start_global, end_global, start_part, end_part;
53  gettimeofday( &start_global, NULL );
54 
55  // Get the whole MPS into left-canonical form
56  gettimeofday( &start_part, NULL );
57  left_normalize( MPS[ L - 2 ], MPS[ L - 1 ] );
58  left_normalize( MPS[ L - 1 ], NULL );
59  gettimeofday( &end_part, NULL );
60  timings[ CHEMPS2_TIME_S_SPLIT ] += ( end_part.tv_sec - start_part.tv_sec ) + 1e-6 * ( end_part.tv_usec - start_part.tv_usec );
61 
62  if ( am_i_master ){
63  if ( do_3rdm ){
64  cout << "****************************************************" << endl;
65  cout << "*** 2-RDM, 3-RDM, and Correlations calculation ***" << endl;
66  cout << "****************************************************" << endl;
67  } else {
68  cout << "********************************************" << endl;
69  cout << "*** 2-RDM and Correlations calculation ***" << endl;
70  cout << "********************************************" << endl;
71  }
72  }
73 
74  // Make the renormalized operators one site further ( one-dot )
75  gettimeofday( &start_part, NULL );
76  updateMovingRightSafe( L - 2 );
77  gettimeofday( &end_part, NULL );
78  timings[ CHEMPS2_TIME_TENS_TOTAL ] += ( end_part.tv_sec - start_part.tv_sec ) + 1e-6 * ( end_part.tv_usec - start_part.tv_usec );
79 
80  // Calculate the 2DM
81  if ( the2DM != NULL ){ delete the2DM; the2DM = NULL; }
82  the2DM = new TwoDM( denBK, Prob );
83 
84  for ( int siteindex = L - 1; siteindex >= 0; siteindex-- ){
85 
86  gettimeofday( &start_part, NULL );
87  // Specific 2-RDM entries are internally added per MPI processes; after which an allreduce is called
88  the2DM->FillSite( MPS[ siteindex ], Ltensors, F0tensors, F1tensors, S0tensors, S1tensors );
89  gettimeofday( &end_part, NULL );
90  timings[ CHEMPS2_TIME_S_SOLVE ] += ( end_part.tv_sec - start_part.tv_sec ) + 1e-6 * ( end_part.tv_usec - start_part.tv_usec );
91 
92  if ( siteindex > 0 ){
93 
94  gettimeofday( &start_part, NULL );
95  right_normalize( MPS[ siteindex - 1 ], MPS[ siteindex ] );
96  gettimeofday( &end_part, NULL );
97  timings[ CHEMPS2_TIME_S_SPLIT ] += ( end_part.tv_sec - start_part.tv_sec ) + 1e-6 * ( end_part.tv_usec - start_part.tv_usec );
98 
99  gettimeofday( &start_part, NULL );
100  updateMovingLeftSafe2DM( siteindex - 1 );
101  gettimeofday( &end_part, NULL );
102  timings[ CHEMPS2_TIME_TENS_TOTAL ] += ( end_part.tv_sec - start_part.tv_sec ) + 1e-6 * ( end_part.tv_usec - start_part.tv_usec );
103  }
104  }
105 
106  #ifdef CHEMPS2_MPI_COMPILATION
107  gettimeofday( &start_part, NULL );
108  the2DM->mpi_allreduce();
109  gettimeofday( &end_part, NULL );
110  timings[ CHEMPS2_TIME_S_SOLVE ] += ( end_part.tv_sec - start_part.tv_sec ) + 1e-6 * ( end_part.tv_usec - start_part.tv_usec );
111  #endif
112 
114 
115  // Trace, energy, and NOON
116  if ( am_i_master ){
117  cout << " N(N-1) = " << denBK->gN() * ( denBK->gN() - 1 ) << endl;
118  cout << " Double trace of DMRG 2-RDM = " << the2DM->trace() << endl;
119  cout << " Econst + 0.5 * trace(2DM-A * Ham) = " << the2DM->energy() << endl;
120  the2DM->print_noon();
121  }
122 
123  // Calculate the 3DM and Correlations
124  if ( the3DM != NULL ){ delete the3DM; the3DM = NULL; }
125  if ( theCorr != NULL ){ delete theCorr; theCorr = NULL; }
126  if ( do_3rdm ){ the3DM = new ThreeDM( denBK, Prob, disk_3rdm ); }
127  theCorr = new Correlations( denBK, Prob, the2DM );
128  if ( am_i_master ){
129  Gtensors = new TensorGYZ*[ L - 1 ];
130  Ytensors = new TensorGYZ*[ L - 1 ];
131  Ztensors = new TensorGYZ*[ L - 1 ];
132  Ktensors = new TensorKM *[ L - 1 ];
133  Mtensors = new TensorKM *[ L - 1 ];
134  }
135  if ( do_3rdm ){
136  tensor_3rdm_a_J0_doublet = new Tensor3RDM****[ L - 1 ];
137  tensor_3rdm_a_J1_doublet = new Tensor3RDM****[ L - 1 ];
138  tensor_3rdm_a_J1_quartet = new Tensor3RDM****[ L - 1 ];
139  tensor_3rdm_b_J0_doublet = new Tensor3RDM****[ L - 1 ];
140  tensor_3rdm_b_J1_doublet = new Tensor3RDM****[ L - 1 ];
141  tensor_3rdm_b_J1_quartet = new Tensor3RDM****[ L - 1 ];
142  tensor_3rdm_c_J0_doublet = new Tensor3RDM****[ L - 1 ];
143  tensor_3rdm_c_J1_doublet = new Tensor3RDM****[ L - 1 ];
144  tensor_3rdm_c_J1_quartet = new Tensor3RDM****[ L - 1 ];
145  tensor_3rdm_d_J0_doublet = new Tensor3RDM****[ L - 1 ];
146  tensor_3rdm_d_J1_doublet = new Tensor3RDM****[ L - 1 ];
147  tensor_3rdm_d_J1_quartet = new Tensor3RDM****[ L - 1 ];
148 
149  // Calculate the leftmost site contribution to the 3-RDM
150  gettimeofday( &start_part, NULL );
151  the3DM->fill_site( MPS[ 0 ], Ltensors, F0tensors, F1tensors, S0tensors, S1tensors,
152  NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL );
153  gettimeofday( &end_part, NULL );
154  timings[ CHEMPS2_TIME_S_SOLVE ] += ( end_part.tv_sec - start_part.tv_sec ) + 1e-6 * ( end_part.tv_usec - start_part.tv_usec );
155  }
156 
157  for ( int siteindex = 1; siteindex < L; siteindex++ ){
158 
159  // Change MPS gauge
160  gettimeofday( &start_part, NULL );
161  left_normalize( MPS[ siteindex - 1 ], MPS[ siteindex ] );
162  gettimeofday( &end_part, NULL );
163  timings[ CHEMPS2_TIME_S_SPLIT ] += ( end_part.tv_sec - start_part.tv_sec ) + 1e-6 * ( end_part.tv_usec - start_part.tv_usec );
164 
165  // Update 2-RDM, 3-RDM, and Correlations tensors
166  gettimeofday( &start_part, NULL );
167  if ( do_3rdm ){ update_safe_3rdm_operators( siteindex ); }
168  updateMovingRightSafe2DM( siteindex - 1 );
169  if ( am_i_master ){ update_correlations_tensors( siteindex ); }
170  gettimeofday( &end_part, NULL );
171  timings[ CHEMPS2_TIME_TENS_TOTAL ] += ( end_part.tv_sec - start_part.tv_sec ) + 1e-6 * ( end_part.tv_usec - start_part.tv_usec );
172 
173  // Calculate Correlation and 3-RDM diagrams. Specific contributions per MPI process. Afterwards an MPI allreduce/bcast is required.
174  gettimeofday( &start_part, NULL );
175  if ( am_i_master ){ theCorr->FillSite( MPS[ siteindex ], Gtensors, Ytensors, Ztensors, Ktensors, Mtensors ); }
176  if ( do_3rdm ){ the3DM->fill_site( MPS[ siteindex ], Ltensors, F0tensors, F1tensors, S0tensors, S1tensors,
177  tensor_3rdm_a_J0_doublet[ siteindex - 1 ], tensor_3rdm_a_J1_doublet[ siteindex - 1 ], tensor_3rdm_a_J1_quartet[ siteindex - 1 ],
178  tensor_3rdm_b_J0_doublet[ siteindex - 1 ], tensor_3rdm_b_J1_doublet[ siteindex - 1 ], tensor_3rdm_b_J1_quartet[ siteindex - 1 ],
179  tensor_3rdm_c_J0_doublet[ siteindex - 1 ], tensor_3rdm_c_J1_doublet[ siteindex - 1 ], tensor_3rdm_c_J1_quartet[ siteindex - 1 ],
180  tensor_3rdm_d_J0_doublet[ siteindex - 1 ], tensor_3rdm_d_J1_doublet[ siteindex - 1 ], tensor_3rdm_d_J1_quartet[ siteindex - 1 ] ); }
181  gettimeofday( &end_part, NULL );
182  timings[ CHEMPS2_TIME_S_SOLVE ] += ( end_part.tv_sec - start_part.tv_sec ) + 1e-6 * ( end_part.tv_usec - start_part.tv_usec );
183  }
184 
185  // Delete the renormalized operators from boundary L-2 and load the ones from boundary L-3
186  gettimeofday( &start_part, NULL );
187  assert( isAllocated[ L - 2 ] == 1 ); // Renormalized operators exist on the last boundary (L-2) and are moving to the right.
188  assert( isAllocated[ L - 3 ] == 0 ); // Renormalized operators do not exist on boundary L-3.
189  deleteTensors( L - 2, true ); isAllocated[ L - 2 ] = 0; // Delete the renormalized operators on the last boundary (L-2).
190  allocateTensors( L - 3, true ); isAllocated[ L - 3 ] = 1; // Create the renormalized operators on boundary L-3.
191  OperatorsOnDisk( L - 3, true, false ); // Load the renormalized operators on boundary L-3.
192  gettimeofday( &end_part, NULL );
193  timings[ CHEMPS2_TIME_TENS_TOTAL ] += ( end_part.tv_sec - start_part.tv_sec ) + 1e-6 * ( end_part.tv_usec - start_part.tv_usec );
194 
195  #ifdef CHEMPS2_MPI_COMPILATION
196  gettimeofday( &start_part, NULL );
197  theCorr->mpi_broadcast();
198  if (do_3rdm){ the3DM->mpi_allreduce(); }
199  gettimeofday( &end_part, NULL );
200  timings[ CHEMPS2_TIME_S_SOLVE ] += ( end_part.tv_sec - start_part.tv_sec ) + 1e-6 * ( end_part.tv_usec - start_part.tv_usec );
201  #endif
202 
203  if ( do_3rdm ){
205  delete_3rdm_operators( L - 1 );
206  delete [] tensor_3rdm_a_J0_doublet;
207  delete [] tensor_3rdm_a_J1_doublet;
208  delete [] tensor_3rdm_a_J1_quartet;
209  delete [] tensor_3rdm_b_J0_doublet;
210  delete [] tensor_3rdm_b_J1_doublet;
211  delete [] tensor_3rdm_b_J1_quartet;
212  delete [] tensor_3rdm_c_J0_doublet;
213  delete [] tensor_3rdm_c_J1_doublet;
214  delete [] tensor_3rdm_c_J1_quartet;
215  delete [] tensor_3rdm_d_J0_doublet;
216  delete [] tensor_3rdm_d_J1_doublet;
217  delete [] tensor_3rdm_d_J1_quartet;
218  }
219 
220  if ( am_i_master ){
221  for ( int previousindex = 0; previousindex < L - 1; previousindex++ ){
222  delete Gtensors[ previousindex ];
223  delete Ytensors[ previousindex ];
224  delete Ztensors[ previousindex ];
225  delete Ktensors[ previousindex ];
226  delete Mtensors[ previousindex ];
227  }
228  delete [] Gtensors;
229  delete [] Ytensors;
230  delete [] Ztensors;
231  delete [] Ktensors;
232  delete [] Mtensors;
233  }
234 
235  gettimeofday( &end_global, NULL );
236  const double elapsed_global = ( end_global.tv_sec - start_global.tv_sec ) + 1e-6 * ( end_global.tv_usec - start_global.tv_usec );
237 
238  if ( am_i_master ){
239  cout << " Single-orbital entropies (Hamiltonian index order is used!) = [ ";
240  for ( int index = 0; index < L - 1; index++ ){ cout << theCorr->SingleOrbitalEntropy_HAM( index ) << " , "; }
241  cout << theCorr->SingleOrbitalEntropy_HAM( L - 1 ) << " ]." << endl;
242  for ( int power = 0; power <= 2; power++ ){
243  cout << " Idistance(" << power << ") = " << theCorr->MutualInformationDistance( (double) power ) << endl;
244  }
245  if ( do_3rdm ){ cout << " N(N-1)(N-2) = " << denBK->gN() * ( denBK->gN() - 1 ) * ( denBK->gN() - 2 ) << endl;
246  cout << " Triple trace of DMRG 3-RDM = " << the3DM->trace() << endl;
247  cout << "***********************************************************" << endl;
248  cout << "*** Timing information 2-RDM, 3-RDM, and Correlations ***" << endl;
249  cout << "***********************************************************" << endl; }
250  else { cout << "***************************************************" << endl;
251  cout << "*** Timing information 2-RDM and Correlations ***" << endl;
252  cout << "***************************************************" << endl; }
253  cout << "*** Elapsed wall time = " << elapsed_global << " seconds" << endl;
254  cout << "*** |--> MPS gauge change = " << timings[ CHEMPS2_TIME_S_SPLIT ] << " seconds" << endl;
255  cout << "*** |--> Diagram calc = " << timings[ CHEMPS2_TIME_S_SOLVE ] << " seconds" << endl;
256  print_tensor_update_performance();
257  if ( do_3rdm ){ cout << "***********************************************************" << endl; }
258  else { cout << "***************************************************" << endl; }
259  }
260 
261 }
262 
263 void CheMPS2::DMRG::print_tensor_update_performance() const{
264 
265  cout << "*** |--> Tensor update = " << timings[ CHEMPS2_TIME_TENS_TOTAL ] << " seconds" << endl;
266  cout << "*** |--> create = " << timings[ CHEMPS2_TIME_TENS_ALLOC ] << " seconds" << endl;
267  cout << "*** |--> destroy = " << timings[ CHEMPS2_TIME_TENS_FREE ] << " seconds" << endl;
268  cout << "*** |--> disk write = " << timings[ CHEMPS2_TIME_DISK_WRITE ] << " seconds" << endl;
269  cout << "*** |--> disk read = " << timings[ CHEMPS2_TIME_DISK_READ ] << " seconds" << endl;
270  cout << "*** |--> calc = " << timings[ CHEMPS2_TIME_TENS_CALC ] << " seconds" << endl;
271  cout << "*** Disk write bandwidth = " << num_double_write_disk * sizeof(double) / ( timings[ CHEMPS2_TIME_DISK_WRITE ] * 1048576 ) << " MB/s" << endl;
272  cout << "*** Disk read bandwidth = " << num_double_read_disk * sizeof(double) / ( timings[ CHEMPS2_TIME_DISK_READ ] * 1048576 ) << " MB/s" << endl;
273 
274 }
275 
276 void CheMPS2::DMRG::left_normalize( TensorT * left_mps, TensorT * right_mps ){
277 
278  #ifdef CHEMPS2_MPI_COMPILATION
279  const bool am_i_master = ( MPIchemps2::mpi_rank() == MPI_CHEMPS2_MASTER );
280  #else
281  const bool am_i_master = true;
282  #endif
283 
284  if ( am_i_master ){
285  const int siteindex = left_mps->gIndex();
286  const SyBookkeeper * theBK = left_mps->gBK();
287  // (J,N,I) = (0,0,0) and (moving_right, prime_last, jw_phase) = (true, true, false)
288  TensorOperator * temp = new TensorOperator( siteindex + 1, 0, 0, 0, true, true, false, theBK, theBK );
289  left_mps->QR( temp );
290  if ( right_mps != NULL ){ right_mps->LeftMultiply( temp ); }
291  delete temp;
292  }
293  #ifdef CHEMPS2_MPI_COMPILATION
294  MPIchemps2::broadcast_tensor( left_mps, MPI_CHEMPS2_MASTER );
295  if ( right_mps != NULL ){ MPIchemps2::broadcast_tensor( right_mps, MPI_CHEMPS2_MASTER ); }
296  #endif
297 
298 }
299 
300 void CheMPS2::DMRG::right_normalize( TensorT * left_mps, TensorT * right_mps ){
301 
302  #ifdef CHEMPS2_MPI_COMPILATION
303  const bool am_i_master = ( MPIchemps2::mpi_rank() == MPI_CHEMPS2_MASTER );
304  #else
305  const bool am_i_master = true;
306  #endif
307 
308  if ( am_i_master ){
309  const int siteindex = right_mps->gIndex();
310  const SyBookkeeper * theBK = right_mps->gBK();
311  // (J,N,I) = (0,0,0) and (moving_right, prime_last, jw_phase) = (true, true, false)
312  TensorOperator * temp = new TensorOperator( siteindex, 0, 0, 0, true, true, false, theBK, theBK );
313  right_mps->LQ( temp );
314  if ( left_mps != NULL ){ left_mps->RightMultiply( temp ); }
315  delete temp;
316  }
317  #ifdef CHEMPS2_MPI_COMPILATION
318  MPIchemps2::broadcast_tensor( right_mps, MPI_CHEMPS2_MASTER );
319  if ( left_mps != NULL ){ MPIchemps2::broadcast_tensor( left_mps, MPI_CHEMPS2_MASTER ); }
320  #endif
321 
322 }
323 
324 double CheMPS2::DMRG::getSpecificCoefficient(int * coeff) const{
325 
326  int * alpha = new int[ L ];
327  int * beta = new int[ L ];
328  for ( int orb=0; orb<L; orb++ ){
329  assert( ( coeff[orb] >= 0 ) && ( coeff[orb] <= 2 ) );
330  if ( coeff[orb] == 0 ){ alpha[orb] = 0; beta[orb] = 0; }
331  if ( coeff[orb] == 1 ){ alpha[orb] = 1; beta[orb] = 0; }
332  if ( coeff[orb] == 2 ){ alpha[orb] = 1; beta[orb] = 1; }
333  }
334  const double FCIcoeff = getFCIcoefficient( alpha, beta );
335  delete [] alpha;
336  delete [] beta;
337  return FCIcoeff;
338 
339 }
340 
341 double CheMPS2::DMRG::getFCIcoefficient(int * alpha, int * beta, const bool mpi_chemps2_master_only) const{
342 
343  //DMRGcoeff = alpha/beta[Hamindex = Prob->gf2(DMRGindex)]
344 
345  //Check if it's possible
346  {
347  int nTot = 0;
348  int twoSz = 0;
349  int iTot = 0;
350  for (int DMRGindex=0; DMRGindex<L; DMRGindex++){
351  const int HamIndex = (Prob->gReorder()) ? Prob->gf2(DMRGindex) : DMRGindex;
352  assert( ( alpha[HamIndex] == 0 ) || ( alpha[HamIndex] == 1 ) );
353  assert( ( beta[HamIndex] == 0 ) || ( beta[HamIndex] == 1 ) );
354  nTot += alpha[HamIndex] + beta[HamIndex];
355  twoSz += alpha[HamIndex] - beta[HamIndex];
356  if ((alpha[HamIndex]+beta[HamIndex])==1){ iTot = Irreps::directProd(iTot,denBK->gIrrep(DMRGindex)); }
357  }
358  if ( Prob->gN() != nTot ){
359  cout << "DMRG::getFCIcoefficient : Ndesired = " << Prob->gN() << " and Ntotal in alpha and beta strings = " << nTot << endl;
360  return 0.0;
361  }
362  // 2Sz can be -Prob->2S() ; -Prob->2S()+2 ; -Prob->2S()+4 ; ... ; Prob->2S()
363  if ( ( Prob->gTwoS() < twoSz ) || ( twoSz < -Prob->gTwoS() ) || ( ( Prob->gTwoS() - twoSz ) % 2 != 0 ) ){
364  cout << "DMRG::getFCIcoefficient : 2Sdesired = " << Prob->gTwoS() << " and 2Sz in alpha and beta strings = " << twoSz << endl;
365  return 0.0;
366  }
367  if ( Prob->gIrrep() != iTot ){
368  cout << "DMRG::getFCIcoefficient : Idesired = " << Prob->gIrrep() << " and Irrep of alpha and beta strings = " << iTot << endl;
369  return 0.0;
370  }
371  }
372 
373  double theCoeff = 2.0; // A FCI coefficient always lies in between -1.0 and 1.0
374  #ifdef CHEMPS2_MPI_COMPILATION
375  if (( MPIchemps2::mpi_rank() == MPI_CHEMPS2_MASTER ) || ( mpi_chemps2_master_only == false ))
376  #endif
377  {
378 
379  //Construct necessary arrays
380  int Dmax = 1;
381  for (int DMRGindex=1; DMRGindex<L; DMRGindex++){
382  const int DtotBound = denBK->gTotDimAtBound(DMRGindex);
383  if (DtotBound>Dmax){ Dmax = DtotBound; }
384  }
385  double * arrayL = new double[Dmax];
386  double * arrayR = new double[Dmax];
387  int * twoSL = new int[L];
388  int * twoSR = new int[L];
389  int * jumpL = new int[L+1];
390  int * jumpR = new int[L+1];
391 
392  //Start the iterator
393  int num_SL = 0;
394  jumpL[num_SL] = 0;
395  int dimFirst = 1;
396  jumpL[num_SL+1] = jumpL[num_SL] + dimFirst;
397  twoSL[num_SL] = 0;
398  num_SL++;
399  arrayL[0] = 1.0;
400  int NL = 0;
401  int IL = 0;
402  int twoSLz = 0;
403 
404  for (int DMRGindex=0; DMRGindex<L; DMRGindex++){
405 
406  //Clear the right array
407  for (int count = 0; count < Dmax; count++){ arrayR[count] = 0.0; }
408 
409  //The local occupation
410  const int HamIndex = (Prob->gReorder()) ? Prob->gf2(DMRGindex) : DMRGindex;
411  const int Nlocal = alpha[HamIndex] + beta[HamIndex];
412  const int twoSzloc = alpha[HamIndex] - beta[HamIndex];
413 
414  //The right symmetry sectors
415  const int NR = NL + Nlocal;
416  const int twoSRz = twoSLz + twoSzloc;
417  const int IR = (( Nlocal == 1 ) ? (Irreps::directProd(IL,denBK->gIrrep(DMRGindex))) : IL);
418 
419  int num_SR = 0;
420  jumpR[num_SR] = 0;
421  const int spread = ( ( Nlocal == 1 ) ? 1 : 0 );
422  for ( int cntSL = 0; cntSL < num_SL; cntSL++ ){
423  for ( int TwoSRattempt = twoSL[cntSL] - spread; TwoSRattempt <= twoSL[cntSL] + spread; TwoSRattempt+=2 ){
424  bool encountered = false;
425  for ( int cntSR = 0; cntSR < num_SR; cntSR++ ){
426  if ( twoSR[cntSR] == TwoSRattempt ){
427  encountered = true;
428  }
429  }
430  if ( encountered == false ){
431  const int dimR = denBK->gCurrentDim(DMRGindex+1,NR,TwoSRattempt,IR);
432  if ( dimR > 0 ){
433  jumpR[num_SR+1] = jumpR[num_SR] + dimR;
434  twoSR[num_SR] = TwoSRattempt;
435  num_SR++;
436  }
437  }
438  }
439  }
440  assert( jumpR[num_SR] <= Dmax );
441 
442  for ( int cntSR = 0; cntSR < num_SR; cntSR++ ){
443  int TwoSRvalue = twoSR[ cntSR ];
444  int dimR = jumpR[ cntSR+1 ] - jumpR[ cntSR ];
445  for ( int TwoSLvalue = TwoSRvalue - spread; TwoSLvalue <= TwoSRvalue + spread; TwoSLvalue += 2 ){
446 
447  int indexSL = -1;
448  for ( int cntSL = 0; cntSL < num_SL; cntSL++ ){
449  if ( twoSL[cntSL] == TwoSLvalue ){
450  indexSL = cntSL;
451  cntSL = num_SL; //exit loop
452  }
453  }
454  if ( indexSL != -1 ){
455  int dimL = jumpL[ indexSL+1 ] - jumpL[ indexSL ];
456  double * Tblock = MPS[DMRGindex]->gStorage(NL,TwoSLvalue,IL,NR,TwoSRvalue,IR);
457  double prefactor = sqrt( TwoSRvalue + 1 )
458  * Wigner::wigner3j(TwoSLvalue, spread, TwoSRvalue, twoSLz, twoSzloc, -twoSRz)
459  * Special::phase( -TwoSLvalue + spread - twoSRz );
460  double add2array = 1.0;
461  char notrans = 'N';
462  dgemm_( &notrans, &notrans, &dimFirst, &dimR, &dimL, &prefactor, arrayL + jumpL[indexSL], &dimFirst, Tblock, &dimL, &add2array, arrayR + jumpR[cntSR], &dimFirst);
463  }
464  }
465  }
466 
467  //Swap L <--> R
468  {
469  double * temp = arrayR;
470  arrayR = arrayL;
471  arrayL = temp;
472  int * temp2 = twoSR;
473  twoSR = twoSL;
474  twoSL = temp2;
475  temp2 = jumpR;
476  jumpR = jumpL;
477  jumpL = temp2;
478  num_SL = num_SR;
479  NL = NR;
480  IL = IR;
481  twoSLz = twoSRz;
482  }
483  }
484 
485  theCoeff = arrayL[0];
486 
487  assert( num_SL == 1 );
488  assert( jumpL[1] == 1 );
489  assert( twoSL[0] == Prob->gTwoS() );
490  assert( NL == Prob->gN() );
491  assert( IL == Prob->gIrrep() );
492 
493  delete [] arrayL;
494  delete [] arrayR;
495  delete [] twoSL;
496  delete [] twoSR;
497  delete [] jumpL;
498  delete [] jumpR;
499 
500  }
501 
502  #ifdef CHEMPS2_MPI_COMPILATION
503  if ( mpi_chemps2_master_only ){ MPIchemps2::broadcast_array_double( &theCoeff, 1, MPI_CHEMPS2_MASTER ); }
504  #endif
505  return theCoeff;
506 
507 }
508 
509 double ** CheMPS2::DMRG::prepare_excitations(Sobject * denS){
510 
511  double ** VeffTilde = new double*[nStates-1];
512  for (int state=0; state<nStates-1; state++){
513  #ifdef CHEMPS2_MPI_COMPILATION
514  if ( MPIchemps2::owner_specific_excitation(L,state) == MPIchemps2::mpi_rank() ){
515  #endif
516  VeffTilde[state] = new double[denS->gKappa2index(denS->gNKappa())];
517  calcVeffTilde(VeffTilde[state], denS, state);
518  #ifdef CHEMPS2_MPI_COMPILATION
519  } else { VeffTilde[state] = NULL; }
520  #endif
521  }
522  return VeffTilde;
523 
524 }
525 
526 void CheMPS2::DMRG::cleanup_excitations(double ** VeffTilde) const{
527 
528  for (int state=0; state<nStates-1; state++){
529  #ifdef CHEMPS2_MPI_COMPILATION
530  if ( MPIchemps2::owner_specific_excitation(L,state) == MPIchemps2::mpi_rank() )
531  #endif
532  {
533  delete [] VeffTilde[state];
534  }
535  }
536  delete [] VeffTilde;
537 
538 }
539 
540 void CheMPS2::DMRG::calcVeffTilde(double * result, Sobject * currentS, int state_number){
541 
542  int dimTot = currentS->gKappa2index(currentS->gNKappa());
543  for (int cnt=0; cnt<dimTot; cnt++){ result[cnt] = 0.0; }
544  int index = currentS->gIndex();
545 
546  const int dimL = std::max(denBK->gMaxDimAtBound(index), Exc_BKs[state_number]->gMaxDimAtBound(index) );
547  const int dimR = std::max(denBK->gMaxDimAtBound(index+2), Exc_BKs[state_number]->gMaxDimAtBound(index+2) );
548  double * workmem = new double[dimL * dimR];
549 
550  //Construct Sup
551  Sobject * Sup = new Sobject( index, Exc_BKs[ state_number ] );
552  Sup->Join(Exc_MPSs[state_number][index],Exc_MPSs[state_number][index+1]);
553 
554  //Construct VeffTilde
555  const double prefactor = sqrt(Exc_Eshifts[state_number]) / (Prob->gTwoS() + 1.0);
556  for (int ikappa=0; ikappa<currentS->gNKappa(); ikappa++){
557  int NL = currentS->gNL(ikappa);
558  int TwoSL = currentS->gTwoSL(ikappa);
559  int IL = currentS->gIL(ikappa);
560  int N1 = currentS->gN1(ikappa);
561  int N2 = currentS->gN2(ikappa);
562  int TwoJ = currentS->gTwoJ(ikappa);
563  int NR = currentS->gNR(ikappa);
564  int TwoSR = currentS->gTwoSR(ikappa);
565  int IR = currentS->gIR(ikappa);
566 
567  //Check if block also exists for other MPS
568  int kappaSup = Sup->gKappa(NL, TwoSL, IL, N1, N2, TwoJ, NR, TwoSR, IR);
569  if (kappaSup!=-1){
570 
571  int dimLdown = denBK->gCurrentDim(index, NL,TwoSL,IL);
572  int dimLup = Exc_BKs[state_number]->gCurrentDim(index, NL,TwoSL,IL);
573  int dimRdown = denBK->gCurrentDim(index+2,NR,TwoSR,IR);
574  int dimRup = Exc_BKs[state_number]->gCurrentDim(index+2,NR,TwoSR,IR);
575 
576  //Do sqrt( (TwoJR+1) * Eshift ) / (TwoStarget+1) times (OL * Sup)_{block} --> workmem
577  double * SupPart = Sup->gStorage() + Sup->gKappa2index(kappaSup);
578  double alpha = prefactor * sqrt(TwoSR+1.0);
579  if (index==0){
580 
581  int dimBlock = dimLup * dimRup;
582  int inc = 1;
583  dcopy_(&dimBlock,SupPart,&inc,workmem,&inc);
584  dscal_(&dimBlock,&alpha,workmem,&inc);
585 
586  } else {
587 
588  char notrans = 'N';
589  double beta = 0.0;
590  double * Opart = Exc_Overlaps[state_number][index-1]->gStorage(NL,TwoSL,IL,NL,TwoSL,IL);
591  dgemm_(&notrans,&notrans,&dimLdown,&dimRup,&dimLup,&alpha,Opart,&dimLdown,SupPart,&dimLup,&beta,workmem,&dimLdown);
592 
593  }
594 
595  //Do (workmem * OR)_{block} --> result + jumpCurrentS
596  int jumpCurrentS = currentS->gKappa2index(ikappa);
597  if (index==L-2){
598 
599  int dimBlock = dimLdown * dimRdown;
600  int inc = 1;
601  dcopy_(&dimBlock, workmem, &inc, result + jumpCurrentS, &inc);
602 
603  } else {
604 
605  char trans = 'T';
606  char notrans = 'N';
607  alpha = 1.0;
608  double beta = 0.0; //set
609  double * Opart = Exc_Overlaps[state_number][index+1]->gStorage(NR,TwoSR,IR,NR,TwoSR,IR);
610  dgemm_(&notrans,&trans,&dimLdown,&dimRdown,&dimRup,&alpha,workmem,&dimLdown,Opart,&dimRdown,&beta,result+jumpCurrentS,&dimLdown);
611 
612  }
613  }
614  }
615 
616  //Deallocate everything
617  delete Sup;
618  delete [] workmem;
619 
620 }
621 
622 void CheMPS2::DMRG::calc_overlaps( const bool moving_right ){
623 
624  for ( int state = 0; state < nStates - 1; state++ ){
625 
626  double overlap;
627  #ifdef CHEMPS2_MPI_COMPILATION
628  const int OWNER = MPIchemps2::owner_specific_excitation( L, state );
629  if ( OWNER == MPIchemps2::mpi_rank() ){
630  #endif
631  if ( moving_right ){
632  TensorO * first = new TensorO( L - 1, moving_right, denBK, Exc_BKs[ state ] );
633  TensorO * second = new TensorO( L, moving_right, denBK, Exc_BKs[ state ] );
634  first->update_ownmem( MPS[ L - 2 ], Exc_MPSs[ state ][ L - 2 ], Exc_Overlaps[ state ][ L - 3 ] );
635  second->update_ownmem( MPS[ L - 1 ], Exc_MPSs[ state ][ L - 1 ], first );
636  overlap = second->gStorage()[ 0 ];
637  delete second;
638  delete first;
639  } else {
640  TensorO * first = new TensorO( 1, moving_right, denBK, Exc_BKs[ state ] );
641  TensorO * second = new TensorO( 0, moving_right, denBK, Exc_BKs[ state ] );
642  first->update_ownmem( MPS[ 1 ], Exc_MPSs[ state ][ 1 ], Exc_Overlaps[ state ][ 1 ] );
643  second->update_ownmem( MPS[ 0 ], Exc_MPSs[ state ][ 0 ], first );
644  overlap = second->gStorage()[ 0 ] / ( Prob->gTwoS() + 1 );
645  delete second;
646  delete first;
647  }
648  #ifdef CHEMPS2_MPI_COMPILATION
649  }
650  MPIchemps2::broadcast_array_double( &overlap, 1, OWNER );
651 
652  if ( MPIchemps2::mpi_rank() == MPI_CHEMPS2_MASTER )
653  #endif
654  { cout << "*** The overlap < " << nStates-1 << " | " << state << " > = " << overlap << endl; }
655 
656  }
657 
658 }
659 
int gTwoSL(const int ikappa) const
Get the left spin symmetry of block ikappa.
Definition: Sobject.cpp:203
void FillSite(TensorT *denT, TensorL ***Ltens, TensorF0 ****F0tens, TensorF1 ****F1tens, TensorS0 ****S0tens, TensorS1 ****S1tens)
Fill the 2DM terms with as second site index denT->gIndex()
Definition: TwoDM.cpp:444
void mpi_allreduce()
Add the 3-RDM elements of all MPI processes.
int gTotDimAtBound(const int boundary) const
Get the total reduced virtual dimension at a certain boundary.
void RightMultiply(Tensor *Mx)
Multiply at the right with a diagonal TensorOperator.
Definition: TensorT.cpp:412
void LQ(Tensor *Lstorage)
Right-normalization.
Definition: TensorT.cpp:291
bool gReorder() const
Get whether the Hamiltonian orbitals are reordered for the DMRG calculation.
Definition: Problem.cpp:346
void correct_higher_multiplicities()
After the whole 3-RDM is filled, a prefactor for higher multiplicities should be applied.
Definition: ThreeDM.cpp:316
int gIndex() const
Get the location index.
Definition: TensorT.cpp:161
void FillSite(TensorT *denT, TensorGYZ **Gtensors, TensorGYZ **Ytensors, TensorGYZ **Ztensors, TensorKM **Ktensors, TensorKM **Mtensors)
Fill at the current step of the iterations the two-orbital mutual information and the remaining part ...
int gCurrentDim(const int boundary, const int N, const int TwoS, const int irrep) const
Get the current virtual dimensions.
int gKappa2index(const int kappa) const
Get the storage jump corresponding to a certain symmetry block.
Definition: Sobject.cpp:190
static int mpi_rank()
Get the rank of this MPI process.
Definition: MPIchemps2.h:131
double SingleOrbitalEntropy_HAM(const int index) const
Get the single-orbital entropy for a certain site, using hte HAM indices.
static int phase(const int power_times_two)
Phase function.
Definition: Special.h:36
void mpi_broadcast()
Broadcast the diradical correlation function and the two-orbital mutual information.
int gTwoJ(const int ikappa) const
Get the central spin symmetry of block ikappa.
Definition: Sobject.cpp:207
double energy() const
Calculate the energy based on the 2DM-A.
Definition: TwoDM.cpp:216
double trace()
Return the trace (should be N(N-1)(N-2))
Definition: ThreeDM.cpp:191
void correct_higher_multiplicities()
After the whole 2-RDM is filled, a prefactor for higher multiplicities should be applied.
Definition: TwoDM.cpp:629
void update_ownmem(TensorT *mps_tensor_up, TensorT *mps_tensor_down, TensorO *previous)
Update the previous TensorO.
Definition: TensorO.cpp:40
int gN() const
Get the targeted particle number.
Definition: Problem.h:68
int gNKappa() const
Get the number of symmetry blocks.
Definition: Sobject.cpp:166
int gTwoS() const
Get twice the targeted spin.
Definition: Problem.h:64
double * gStorage()
Get the pointer to the storage.
double getFCIcoefficient(int *alpha, int *beta, const bool mpi_chemps2_master_only=true) const
Get a specific FCI coefficient. The arrays alpha and beta contain the alpha and beta occupation numbe...
int gIL(const int ikappa) const
Get the left irrep symmetry of block ikappa.
Definition: Sobject.cpp:204
double * gStorage()
Get the pointer to the storage.
Definition: Sobject.cpp:168
int gIrrep(const int nOrb) const
Get an orbital irrep.
Definition: Problem.cpp:336
double trace() const
Return the double trace of 2DM-A (should be N(N-1))
Definition: TwoDM.cpp:204
double * gStorage()
Get the pointer to the storage.
Definition: TensorT.cpp:134
static double wigner3j(const int two_ja, const int two_jb, const int two_jc, const int two_ma, const int two_mb, const int two_mc)
Wigner-3j symbol (gsl api)
Definition: Wigner.cpp:243
double MutualInformationDistance(const double power) const
Return Idistance(power) (see return for definition)
int gKappa(const int NL, const int TwoSL, const int IL, const int N1, const int N2, const int TwoJ, const int NR, const int TwoSR, const int IR) const
Get the index corresponding to a certain symmetry block.
Definition: Sobject.cpp:172
void Join(TensorT *Tleft, TensorT *Tright)
Join two sites to form a composite S-object.
Definition: Sobject.cpp:212
static int directProd(const int Irrep1, const int Irrep2)
Get the direct product of the irreps with numbers Irrep1 and Irrep2: a bitwise XOR for psi4&#39;s convent...
Definition: Irreps.h:123
void QR(Tensor *Rstorage)
Left-normalization.
Definition: TensorT.cpp:188
void LeftMultiply(Tensor *Mx)
Multiply at the left with a diagonal TensorOperator.
Definition: TensorT.cpp:391
const SyBookkeeper * gBK() const
Get the pointer to the symmetry bookkeeper.
Definition: TensorT.cpp:163
void mpi_allreduce()
Add the 2-RDM elements of all MPI processes.
double getSpecificCoefficient(int *coeff) const
Get a specific FCI coefficient. The array coeff contains the occupation numbers of the L Hamiltonian ...
void print_noon() const
Print the natural orbital occupation numbers.
Definition: TwoDM.cpp:233
int gNL(const int ikappa) const
Get the left particle number symmetry of block ikappa.
Definition: Sobject.cpp:202
void fill_site(TensorT *denT, TensorL ***Ltens, TensorF0 ****F0tens, TensorF1 ****F1tens, TensorS0 ****S0tens, TensorS1 ****S1tens, Tensor3RDM ****dm3_a_J0_doublet, Tensor3RDM ****dm3_a_J1_doublet, Tensor3RDM ****dm3_a_J1_quartet, Tensor3RDM ****dm3_b_J0_doublet, Tensor3RDM ****dm3_b_J1_doublet, Tensor3RDM ****dm3_b_J1_quartet, Tensor3RDM ****dm3_c_J0_doublet, Tensor3RDM ****dm3_c_J1_doublet, Tensor3RDM ****dm3_c_J1_quartet, Tensor3RDM ****dm3_d_J0_doublet, Tensor3RDM ****dm3_d_J1_doublet, Tensor3RDM ****dm3_d_J1_quartet)
Fill the 3-RDM terms corresponding to site denT->gIndex()
Definition: ThreeDM.cpp:398
int gIR(const int ikappa) const
Get the right irrep symmetry of block ikappa.
Definition: Sobject.cpp:210
int gTwoSR(const int ikappa) const
Get the right spin symmetry of block ikappa.
Definition: Sobject.cpp:209
int gIndex() const
Get the location index.
Definition: Sobject.cpp:200
int gMaxDimAtBound(const int boundary) const
Get the maximum virtual dimension at a certain boundary.
void calc_rdms_and_correlations(const bool do_3rdm, const bool disk_3rdm=false)
Calculate the reduced density matrices and correlations. Afterwards the MPS is again in LLLLLLLC gaug...
int gNR(const int ikappa) const
Get the right particle number symmetry of block ikappa.
Definition: Sobject.cpp:208
int gIrrep(const int orbital) const
Get an orbital irrep.
int gN1(const int ikappa) const
Get the left local particle number symmetry of block ikappa.
Definition: Sobject.cpp:205
int gN() const
Get the targeted particle number.
int gf2(const int DMRGOrb) const
Get the Ham index corresponding to a DMRG index.
Definition: Problem.cpp:348
int gN2(const int ikappa) const
Get the right local particle number symmetry of block ikappa.
Definition: Sobject.cpp:206