CheMPS2
Correlations.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 <stdlib.h>
21 #include <iostream>
22 #include <sstream>
23 #include <algorithm>
24 #include <math.h>
25 
26 #include "Correlations.h"
27 #include "Lapack.h"
28 #include "MPIchemps2.h"
29 
30 using std::cout;
31 using std::endl;
32 using std::max;
33 using std::min;
34 
35 CheMPS2::Correlations::Correlations(const SyBookkeeper * denBKIn, const Problem * ProbIn, TwoDM * the2DMin){
36 
37  denBK = denBKIn;
38  Prob = ProbIn;
39  the2DM = the2DMin;
40 
41  L = denBK->gL();
42 
43  Cspin = new double[L*L];
44  Cdens = new double[L*L];
45  Cspinflip = new double[L*L];
46  Cdirad = new double[L*L];
47  MutInfo = new double[L*L];
48 
49  for (int cnt=0; cnt<L*L; cnt++){ Cspin[cnt] = 0.0; }
50  for (int cnt=0; cnt<L*L; cnt++){ Cdens[cnt] = 0.0; }
51  for (int cnt=0; cnt<L*L; cnt++){ Cspinflip[cnt] = 0.0; }
52  for (int cnt=0; cnt<L*L; cnt++){ Cdirad[cnt] = 0.0; }
53  for (int cnt=0; cnt<L*L; cnt++){ MutInfo[cnt] = 0.0; }
54 
55  FillSpinDensSpinflip();
56 
57 }
58 
60 
61  delete [] Cspin;
62  delete [] Cdens;
63  delete [] Cspinflip;
64  delete [] Cdirad;
65  delete [] MutInfo;
66 
67 }
68 
69 void CheMPS2::Correlations::FillSpinDensSpinflip(){
70 
71  //Spin
72  for (int row=0; row<L; row++){
73  for (int col=0; col<L; col++){
74  Cspin[row + L*col] = the2DM->getTwoDMB_DMRG(row,col,row,col);
75  }
76  Cspin[row + L*row] += the2DM->get1RDM_DMRG(row, row);
77  }
78 
79  //Density
80  for (int row=0; row<L; row++){
81  for (int col=0; col<L; col++){
82  Cdens[row + L*col] = the2DM->getTwoDMA_DMRG(row,col,row,col) - the2DM->get1RDM_DMRG(row, row) * the2DM->get1RDM_DMRG(col, col);
83  }
84  Cdens[row + L*row] += the2DM->get1RDM_DMRG(row, row);
85  }
86 
87  //Spin-flip
88  for (int row=0; row<L; row++){
89  for (int col=0; col<L; col++){
90  Cspinflip[row + L*col] = 0.5 * ( the2DM->getTwoDMB_DMRG(row,col,col,row) - the2DM->getTwoDMA_DMRG(row,col,col,row) );
91  }
92  Cspinflip[row + L*row] += the2DM->get1RDM_DMRG(row, row);
93  }
94 
95  //Singlet diradical: partial fill
96  for (int row=0; row<L; row++){
97  for (int col=0; col<L; col++){
98  Cdirad[row + L*col] = - 0.5 * ( the2DM->get1RDM_DMRG(row, row) - the2DM->getTwoDMA_DMRG(row,row,row,row) )
99  * ( the2DM->get1RDM_DMRG(col, col) - the2DM->getTwoDMA_DMRG(col,col,col,col) );
100  }
101  }
102 
103 }
104 
105 double CheMPS2::Correlations::getCspin_DMRG(const int row, const int col) const{ return Cspin[row + L*col]; }
106 
107 double CheMPS2::Correlations::getCspin_HAM(const int row, const int col) const{
108 
109  //Prob assumes you use DMRG orbs... f1 converts HAM orbs to DMRG orbs
110  if ( Prob->gReorder() ){
111  return getCspin_DMRG( Prob->gf1(row), Prob->gf1(col) );
112  }
113  return getCspin_DMRG( row, col );
114 
115 }
116 
117 double CheMPS2::Correlations::getCdens_DMRG(const int row, const int col) const{ return Cdens[row + L*col]; }
118 
119 double CheMPS2::Correlations::getCdens_HAM(const int row, const int col) const{
120 
121  //Prob assumes you use DMRG orbs... f1 converts HAM orbs to DMRG orbs
122  if ( Prob->gReorder() ){
123  return getCdens_DMRG( Prob->gf1(row), Prob->gf1(col) );
124  }
125  return getCdens_DMRG( row, col );
126 
127 }
128 
129 double CheMPS2::Correlations::getCspinflip_DMRG(const int row, const int col) const{ return Cspinflip[row + L*col]; }
130 
131 double CheMPS2::Correlations::getCspinflip_HAM(const int row, const int col) const{
132 
133  //Prob assumes you use DMRG orbs... f1 converts HAM orbs to DMRG orbs
134  if ( Prob->gReorder() ){
135  return getCspinflip_DMRG( Prob->gf1(row), Prob->gf1(col) );
136  }
137  return getCspinflip_DMRG( row, col );
138 
139 }
140 
141 double CheMPS2::Correlations::getCdirad_DMRG(const int row, const int col) const{ return Cdirad[row + L*col]; }
142 
143 double CheMPS2::Correlations::getCdirad_HAM(const int row, const int col) const{
144 
145  //Prob assumes you use DMRG orbs... f1 converts HAM orbs to DMRG orbs
146  if ( Prob->gReorder() ){
147  return getCdirad_DMRG( Prob->gf1(row), Prob->gf1(col) );
148  }
149  return getCdirad_DMRG( row, col );
150 
151 }
152 
153 double CheMPS2::Correlations::getMutualInformation_DMRG(const int row, const int col) const{ return MutInfo[row + L*col]; }
154 
155 double CheMPS2::Correlations::getMutualInformation_HAM(const int row, const int col) const{
156 
157  //Prob assumes you use DMRG orbs... f1 converts HAM orbs to DMRG orbs
158  if ( Prob->gReorder() ){
159  return getMutualInformation_DMRG( Prob->gf1(row), Prob->gf1(col) );
160  }
161  return getMutualInformation_DMRG( row, col );
162 
163 }
164 
166 
167  const double val4 = 0.5 * the2DM->getTwoDMA_DMRG(index,index,index,index);
168  const double val23 = 0.5 * ( the2DM->get1RDM_DMRG(index, index) - the2DM->getTwoDMA_DMRG(index,index,index,index) );
169  const double val1 = 1.0 - val4 - 2*val23;
170 
171  double entropy = 0.0;
172  if (val1 > CheMPS2::CORRELATIONS_discardEig){ entropy -= val1 * log(val1 ); }
173  if (val23 > CheMPS2::CORRELATIONS_discardEig){ entropy -= 2 * val23 * log(val23); }
174  if (val4 > CheMPS2::CORRELATIONS_discardEig){ entropy -= val4 * log(val4 ); }
175  return entropy;
176 
177 }
178 
180 
181  if ( Prob->gReorder() ){
182  return SingleOrbitalEntropy_DMRG( Prob->gf1(index) );
183  }
184  return SingleOrbitalEntropy_DMRG( index );
185 
186 }
187 
188 double CheMPS2::Correlations::MutualInformationDistance(const double power) const{
189 
190  double Idist = 0.0;
191  for (int row=0; row<L; row++){
192  for (int col=0; col<L; col++){
193  if (row != col){
194  Idist += MutInfo[row + L*col] * pow(abs(row - col), power);
195  }
196  }
197  }
198  return Idist;
199 
200 }
201 
202 #ifdef CHEMPS2_MPI_COMPILATION
204 
205  int size = L*L;
206  MPIchemps2::broadcast_array_double( MutInfo, size, MPI_CHEMPS2_MASTER );
207  MPIchemps2::broadcast_array_double( Cdirad, size, MPI_CHEMPS2_MASTER );
208 
209 }
210 #endif
211 
212 void CheMPS2::Correlations::FillSite(TensorT * denT, TensorGYZ ** Gtensors, TensorGYZ ** Ytensors, TensorGYZ ** Ztensors, TensorKM ** Ktensors, TensorKM ** Mtensors){
213 
214  const int theindex = denT->gIndex();
215  const int MAXDIM = max(denBK->gMaxDimAtBound(theindex), denBK->gMaxDimAtBound(theindex+1));
216  double * workmem = new double[MAXDIM*MAXDIM];
217  int lindimRDM = 16;
218  double * RDM = new double[lindimRDM*lindimRDM];
219  int lwork = 3*lindimRDM - 1;
220  double * work = new double[lwork];
221  double * eigs = new double[lindimRDM];
222 
223  const double prefactorSpin = 1.0/(Prob->gTwoS() + 1.0);
224  const double sqrt_one_half = sqrt(0.5);
225 
226  for (int previousindex=0; previousindex<theindex; previousindex++){
227 
228  const bool equalIrreps = (denBK->gIrrep(previousindex) == denBK->gIrrep(theindex)) ? true : false;
229  const double diag1 = diagram3(denT, Gtensors[previousindex], workmem) * prefactorSpin * 0.5 * sqrt_one_half;
230  const double diag2 = 0.125 * ( the2DM->getTwoDMB_DMRG(previousindex,theindex,theindex,previousindex)
231  - the2DM->getTwoDMA_DMRG(previousindex,theindex,theindex,previousindex) );
232 
233  const double val1 = diagram1(denT, Ytensors[previousindex], workmem) * prefactorSpin; //1x1 block N=0, Sz=0
234  const double val2 = diagram2(denT, Ztensors[previousindex], workmem) * prefactorSpin; //1x1 block N=4, Sz=0
235  const double val3 = diag1 + diag2; //1x1 block N=2, Sz=2*sigma
236  const double val4 = diagram1(denT, Gtensors[previousindex], workmem) * prefactorSpin * sqrt_one_half; //2x2 block N=1, alpha_LL
237  const double val5 = diagram3(denT, Ytensors[previousindex], workmem) * prefactorSpin * 0.5; //2x2 block N=1, alpha_RR
238  const double val6 = (equalIrreps) ? ( diagram4(denT, Ktensors[previousindex], workmem) * prefactorSpin * 0.5 ) : 0.0 ; //2x2 block N=1, alpha_LR
239  const double val7 = diagram2(denT, Gtensors[previousindex], workmem) * prefactorSpin * sqrt_one_half; //2x2 block N=3, alpha_LL
240  const double val8 = diagram3(denT, Ztensors[previousindex], workmem) * prefactorSpin * 0.5; //2x2 block N=3, alpha_RR
241  const double val9 = (equalIrreps) ? ( diagram5(denT, Mtensors[previousindex], workmem) * prefactorSpin * 0.5 ) : 0.0 ; //2x2 block N=3, alpha_LR
242 
243  //4x4 block N=2, Sz=0
244  const double alpha = diagram2(denT, Ytensors[previousindex], workmem) * prefactorSpin;
245  const double gamma = diagram1(denT, Ztensors[previousindex], workmem) * prefactorSpin;
246  const double beta = diag1 - diag2;
247  const double lambda = 2*diag2;
248  const double delta = (equalIrreps) ? ( - diagram5(denT, Ktensors[previousindex], workmem) * prefactorSpin * 0.5 ) : 0.0;
249  const double epsilon = (equalIrreps) ? ( diagram4(denT, Mtensors[previousindex], workmem) * prefactorSpin * 0.5 ) : 0.0;
250  const double kappa = 0.5 * the2DM->getTwoDMA_DMRG(previousindex,previousindex,theindex,theindex);
251 
252  /*
253 
254  [ val1 ]
255  [ val4 val6 ]
256  [ val6 val5 ]
257  [ val4 val6 ]
258  [ val6 val5 ]
259  [ val3 ]
260  [ alpha delta -delta kappa ]
261  [ delta beta lambda epsilon ]
262  [ -delta lambda beta -epsilon ]
263  [ kappa epsilon -epsilon gamma ]
264  [ val3 ]
265  [ val7 val9 ]
266  [ val9 val8 ]
267  [ val7 val9 ]
268  [ val9 val8 ]
269  [ val2 ]
270 
271  */
272 
273  for (int cnt=0; cnt<lindimRDM*lindimRDM; cnt++){ RDM[cnt] = 0.0; }
274  RDM[0 + lindimRDM * 0 ] = val1;
275  RDM[15 + lindimRDM * 15] = val2;
276  RDM[5 + lindimRDM * 5 ] = RDM[10 + lindimRDM * 10] = val3;
277  RDM[1 + lindimRDM * 1 ] = RDM[3 + lindimRDM * 3 ] = val4;
278  RDM[2 + lindimRDM * 2 ] = RDM[4 + lindimRDM * 4 ] = val5;
279  RDM[1 + lindimRDM * 2 ] = RDM[2 + lindimRDM * 1 ] = RDM[3 + lindimRDM * 4 ] = RDM[4 + lindimRDM * 3 ] = val6;
280  RDM[11 + lindimRDM * 11] = RDM[13 + lindimRDM * 13] = val7;
281  RDM[12 + lindimRDM * 12] = RDM[14 + lindimRDM * 14] = val8;
282  RDM[11 + lindimRDM * 12] = RDM[12 + lindimRDM * 11] = RDM[13 + lindimRDM*14] = RDM[14 + lindimRDM*13] = val9;
283  RDM[6 + lindimRDM * 6 ] = alpha;
284  RDM[7 + lindimRDM * 7 ] = RDM[8 + lindimRDM * 8 ] = beta;
285  RDM[9 + lindimRDM * 9 ] = gamma;
286  RDM[6 + lindimRDM * 7 ] = RDM[7 + lindimRDM * 6 ] = delta;
287  RDM[6 + lindimRDM * 8 ] = RDM[8 + lindimRDM * 6 ] = - delta;
288  RDM[7 + lindimRDM * 9 ] = RDM[9 + lindimRDM * 7 ] = epsilon;
289  RDM[8 + lindimRDM * 9 ] = RDM[9 + lindimRDM * 8 ] = - epsilon;
290  RDM[6 + lindimRDM * 9 ] = RDM[9 + lindimRDM * 6 ] = kappa;
291  RDM[7 + lindimRDM * 8 ] = RDM[8 + lindimRDM * 7 ] = lambda;
292 
293  if (CheMPS2::CORRELATIONS_debugPrint){
294  const double RDM_1orb_prev_4 = 0.5 * the2DM->getTwoDMA_DMRG(previousindex,previousindex,previousindex,previousindex);
295  const double RDM_1orb_prev_23 = 0.5 * ( the2DM->get1RDM_DMRG(previousindex, previousindex)
296  - the2DM->getTwoDMA_DMRG(previousindex,previousindex,previousindex,previousindex) );
297  const double RDM_1orb_prev_1 = 1.0 - RDM_1orb_prev_4 - 2*RDM_1orb_prev_23;
298 
299  const double RDM_1orb_curr_4 = 0.5 * the2DM->getTwoDMA_DMRG(theindex,theindex,theindex,theindex);
300  const double RDM_1orb_curr_23 = 0.5 * ( the2DM->get1RDM_DMRG(theindex, theindex) - the2DM->getTwoDMA_DMRG(theindex,theindex,theindex,theindex) );
301  const double RDM_1orb_curr_1 = 1.0 - RDM_1orb_curr_4 - 2*RDM_1orb_curr_23;
302 
303  cout << " Correlations::FillSite : Looking at DMRG sites (" << previousindex << "," << theindex << ")." << endl;
304 
305  //Check 1 : full trace
306  double fulltrace = 0.0;
307  for (int cnt=0; cnt<lindimRDM; cnt++){ fulltrace += RDM[ cnt * ( 1 + lindimRDM ) ]; }
308  cout << " 1 - trace(2-orb RDM) = " << 1.0 - fulltrace << endl;
309 
310  //Check 2 : trace over orb previousindex
311  double getal1 = RDM[0 + lindimRDM * 0 ] + RDM[1 + lindimRDM * 1 ] + RDM[3 + lindimRDM * 3 ] + RDM[9 + lindimRDM * 9 ] - RDM_1orb_curr_1;
312  double getal2 = RDM[2 + lindimRDM * 2 ] + RDM[5 + lindimRDM * 5 ] + RDM[8 + lindimRDM * 8 ] + RDM[12 + lindimRDM * 12] - RDM_1orb_curr_23;
313  double getal3 = RDM[4 + lindimRDM * 4 ] + RDM[7 + lindimRDM * 7 ] + RDM[10 + lindimRDM * 10] + RDM[14 + lindimRDM * 14] - RDM_1orb_curr_23;
314  double getal4 = RDM[6 + lindimRDM * 6 ] + RDM[11 + lindimRDM * 11] + RDM[13 + lindimRDM * 13] + RDM[15 + lindimRDM * 15] - RDM_1orb_curr_4;
315  double RMS = sqrt(getal1*getal1 + getal2*getal2 + getal3*getal3 + getal4*getal4);
316  cout << " 2-norm difference of one-orb RDM of the CURRENT site via 2DM and via trace(2-orb RDM) = " << RMS << endl;
317 
318  //Check 3 : trace over orb currentindex
319  getal1 = RDM[0 + lindimRDM * 0 ] + RDM[2 + lindimRDM * 2 ] + RDM[4 + lindimRDM * 4 ] + RDM[6 + lindimRDM * 6 ] - RDM_1orb_prev_1;
320  getal2 = RDM[1 + lindimRDM * 1 ] + RDM[5 + lindimRDM * 5 ] + RDM[7 + lindimRDM * 7 ] + RDM[11 + lindimRDM * 11] - RDM_1orb_prev_23;
321  getal3 = RDM[3 + lindimRDM * 3 ] + RDM[8 + lindimRDM * 8 ] + RDM[10 + lindimRDM * 10] + RDM[13 + lindimRDM * 13] - RDM_1orb_prev_23;
322  getal4 = RDM[9 + lindimRDM * 9 ] + RDM[12 + lindimRDM * 12] + RDM[14 + lindimRDM * 14] + RDM[15 + lindimRDM * 15] - RDM_1orb_prev_4;
323  RMS = sqrt(getal1*getal1 + getal2*getal2 + getal3*getal3 + getal4*getal4);
324  cout << " 2-norm difference of one-orb RDM of the PREVIOUS site via 2DM and via trace(2-orb RDM) = " << RMS << endl;
325  }
326 
327  char jobz = 'N'; //eigenvalues only
328  char uplo = 'U';
329  int info;
330  dsyev_(&jobz, &uplo, &lindimRDM, RDM, &lindimRDM, eigs, work, &lwork, &info);
331 
332  double entropy = 0.0;
333  for (int cnt=0; cnt<lindimRDM; cnt++){
334  if (eigs[cnt] > CheMPS2::CORRELATIONS_discardEig){ //With discardEig = 1e-100, the discarded contributions are smaller than 2.4e-98
335  entropy -= eigs[cnt] * log(eigs[cnt]);
336  }
337  }
338 
339  const double thisMutInfo = 0.5 * ( SingleOrbitalEntropy_DMRG(previousindex) + SingleOrbitalEntropy_DMRG(theindex) - entropy );
340  MutInfo[previousindex + L * theindex] = MutInfo[theindex + L * previousindex] = thisMutInfo;
341  Cdirad[previousindex + L * theindex] += 2 * beta;
342  Cdirad[theindex + L * previousindex] += 2 * beta;
343 
344  }
345 
346  delete [] eigs;
347  delete [] work;
348  delete [] RDM;
349  delete [] workmem;
350 
351 }
352 
353 double CheMPS2::Correlations::diagram1(TensorT * denT, TensorGYZ * denY, double * workmem) const{ //checked
354 
355  const int theindex = denT->gIndex();
356 
357  double total = 0.0;
358 
359  for (int NR = denBK->gNmin(theindex+1); NR <= denBK->gNmax(theindex+1); NR++){
360  for (int TwoSR = denBK->gTwoSmin(theindex+1,NR); TwoSR <= denBK->gTwoSmax(theindex+1,NR); TwoSR += 2){
361  for (int IR = 0; IR < denBK->getNumberOfIrreps(); IR++){
362 
363  int dimL = denBK->gCurrentDim(theindex, NR, TwoSR, IR);
364  int dimR = denBK->gCurrentDim(theindex+1, NR, TwoSR, IR);
365 
366  if ((dimL>0)&&(dimR>0)){
367 
368  double * blockT = denT->gStorage(NR,TwoSR,IR,NR,TwoSR,IR);
369  double * blockY = denY->gStorage(NR,TwoSR,IR,NR,TwoSR,IR);
370 
371  char notrans = 'N';
372  double alpha = 1.0;
373  double beta = 0.0; //SET
374  dgemm_(&notrans, &notrans, &dimL, &dimR, &dimL, &alpha, blockY, &dimL, blockT, &dimL, &beta, workmem, &dimL);
375 
376  int length = dimL*dimR;
377  int inc = 1;
378  total += (TwoSR + 1.0) * ddot_(&length, workmem, &inc, blockT, &inc);
379 
380  }
381  }
382  }
383  }
384 
385  return total;
386 
387 }
388 
389 double CheMPS2::Correlations::diagram2(TensorT * denT, TensorGYZ * denZ, double * workmem) const{ //checked
390 
391  const int theindex = denT->gIndex();
392 
393  double total = 0.0;
394 
395  for (int NR = denBK->gNmin(theindex+1); NR <= denBK->gNmax(theindex+1); NR++){
396  for (int TwoSR = denBK->gTwoSmin(theindex+1,NR); TwoSR <= denBK->gTwoSmax(theindex+1,NR); TwoSR += 2){
397  for (int IR = 0; IR < denBK->getNumberOfIrreps(); IR++){
398 
399  int dimL = denBK->gCurrentDim(theindex, NR-2, TwoSR, IR);
400  int dimR = denBK->gCurrentDim(theindex+1, NR, TwoSR, IR);
401 
402  if ((dimL>0)&&(dimR>0)){
403 
404  double * blockT = denT->gStorage(NR-2,TwoSR,IR,NR, TwoSR,IR);
405  double * blockZ = denZ->gStorage(NR-2,TwoSR,IR,NR-2,TwoSR,IR);
406 
407  char notrans = 'N';
408  double alpha = 1.0;
409  double beta = 0.0; //SET
410  dgemm_(&notrans, &notrans, &dimL, &dimR, &dimL, &alpha, blockZ, &dimL, blockT, &dimL, &beta, workmem, &dimL);
411 
412  int length = dimL*dimR;
413  int inc = 1;
414  total += (TwoSR + 1.0) * ddot_(&length, workmem, &inc, blockT, &inc);
415 
416  }
417  }
418  }
419  }
420 
421  return total;
422 
423 }
424 
425 double CheMPS2::Correlations::diagram3(TensorT * denT, TensorGYZ * denG, double * workmem) const{ //checked
426 
427  const int theindex = denT->gIndex();
428 
429  double total = 0.0;
430 
431  for (int NR = denBK->gNmin(theindex+1); NR <= denBK->gNmax(theindex+1); NR++){
432  for (int TwoSR = denBK->gTwoSmin(theindex+1,NR); TwoSR <= denBK->gTwoSmax(theindex+1,NR); TwoSR += 2){
433  for (int IR = 0; IR < denBK->getNumberOfIrreps(); IR++){
434 
435  int dimR = denBK->gCurrentDim(theindex+1, NR, TwoSR, IR);
436  const int IL = Irreps::directProd( IR, denBK->gIrrep(theindex) );
437 
438  if (dimR>0){
439 
440  for (int TwoSL = TwoSR-1; TwoSL <= TwoSR+1; TwoSL+=2){
441 
442  int dimL = denBK->gCurrentDim(theindex, NR-1, TwoSL, IL);
443 
444  if (dimL>0){
445 
446  double * blockT = denT->gStorage(NR-1,TwoSL,IL,NR, TwoSR,IR);
447  double * blockG = denG->gStorage(NR-1,TwoSL,IL,NR-1,TwoSL,IL);
448 
449  char notrans = 'N';
450  double alpha = 1.0;
451  double beta = 0.0; //SET
452  dgemm_(&notrans, &notrans, &dimL, &dimR, &dimL, &alpha, blockG, &dimL, blockT, &dimL, &beta, workmem, &dimL);
453 
454  int length = dimL*dimR;
455  int inc = 1;
456  total += (TwoSR + 1.0) * ddot_(&length, workmem, &inc, blockT, &inc);
457 
458  }
459  }
460  }
461  }
462  }
463  }
464 
465  return total;
466 
467 }
468 
469 double CheMPS2::Correlations::diagram4(TensorT * denT, TensorKM * denK, double * workmem) const{ //checked
470 
471  const int theindex = denT->gIndex();
472 
473  double total = 0.0;
474 
475  for (int NR = denBK->gNmin(theindex+1); NR <= denBK->gNmax(theindex+1); NR++){
476  for (int TwoSR = denBK->gTwoSmin(theindex+1,NR); TwoSR <= denBK->gTwoSmax(theindex+1,NR); TwoSR += 2){
477  for (int IR = 0; IR < denBK->getNumberOfIrreps(); IR++){
478 
479  int dimR = denBK->gCurrentDim(theindex+1, NR, TwoSR, IR);
480  int dimLup = denBK->gCurrentDim(theindex, NR, TwoSR, IR);
481  const int ILdown = Irreps::directProd( IR, denBK->gIrrep(theindex) );
482 
483  if ((dimR>0) && (dimLup>0)){
484 
485  for (int TwoSLdown = TwoSR-1; TwoSLdown <= TwoSR+1; TwoSLdown += 2){
486 
487  int dimLdown = denBK->gCurrentDim(theindex, NR-1, TwoSLdown, ILdown);
488 
489  if (dimLdown>0){
490 
491  double * blockTup = denT->gStorage(NR, TwoSR, IR, NR, TwoSR, IR);
492  double * blockTdown = denT->gStorage(NR-1, TwoSLdown, ILdown, NR, TwoSR, IR);
493  double * blockK = denK->gStorage(NR-1, TwoSLdown, ILdown, NR, TwoSR, IR);
494 
495  char notrans = 'N';
496  double alpha = 1.0;
497  double beta = 0.0; //SET
498  dgemm_(&notrans, &notrans, &dimLdown, &dimR, &dimLup, &alpha, blockK, &dimLdown, blockTup, &dimLup, &beta, workmem, &dimLdown);
499 
500  int length = dimLdown*dimR;
501  int inc = 1;
502  total += (TwoSR + 1.0) * ddot_(&length, workmem, &inc, blockTdown, &inc);
503 
504  }
505  }
506  }
507  }
508  }
509  }
510 
511  return total;
512 
513 }
514 
515 double CheMPS2::Correlations::diagram5(TensorT * denT, TensorKM * denM, double * workmem) const{ //checked
516 
517  const int theindex = denT->gIndex();
518 
519  double total = 0.0;
520 
521  for (int NR = denBK->gNmin(theindex+1); NR <= denBK->gNmax(theindex+1); NR++){
522  for (int TwoSR = denBK->gTwoSmin(theindex+1,NR); TwoSR <= denBK->gTwoSmax(theindex+1,NR); TwoSR += 2){
523  for (int IR = 0; IR < denBK->getNumberOfIrreps(); IR++){
524 
525  int dimR = denBK->gCurrentDim(theindex+1, NR, TwoSR, IR);
526  int dimLdown = denBK->gCurrentDim(theindex, NR-2, TwoSR, IR);
527  const int ILup = Irreps::directProd( IR, denBK->gIrrep(theindex) );
528 
529  if ((dimR>0) && (dimLdown>0)){
530 
531  for (int TwoSLup = TwoSR-1; TwoSLup <= TwoSR+1; TwoSLup += 2){
532 
533  int dimLup = denBK->gCurrentDim(theindex, NR-1, TwoSLup, ILup);
534 
535  if (dimLup>0){
536 
537  double * blockTdown = denT->gStorage(NR-2, TwoSR, IR, NR, TwoSR, IR);
538  double * blockTup = denT->gStorage(NR-1, TwoSLup, ILup, NR, TwoSR, IR);
539  double * blockM = denM->gStorage(NR-2, TwoSR, IR, NR-1, TwoSLup, ILup);
540 
541  char notrans = 'N';
542  double alpha = 1.0;
543  double beta = 0.0; //SET
544  dgemm_(&notrans, &notrans, &dimLdown, &dimR, &dimLup, &alpha, blockM, &dimLdown, blockTup, &dimLup, &beta, workmem, &dimLdown);
545 
546  int length = dimLdown*dimR;
547  int inc = 1;
548  const int fase = ((((TwoSLup + 1 - TwoSR)/2)%2)!=0)?-1:1;
549  total += fase * sqrt((TwoSLup+1.0)*(TwoSR+1)) * ddot_(&length, workmem, &inc, blockTdown, &inc);
550 
551  }
552  }
553  }
554  }
555  }
556  }
557 
558  return total;
559 
560 }
561 
562 void CheMPS2::Correlations::PrintTableNice(const double * table, const int sPrecision, const int columnsPerLine) const{
563 
564  std::stringstream thestream;
565  thestream.precision(sPrecision);
566  thestream << std::fixed;
567 
568  int numGroups = L / columnsPerLine;
569  if (numGroups * columnsPerLine < L){ numGroups++; }
570 
571  std::string prefix = " ";
572 
573  for (int groups=0; groups<numGroups; groups++){
574  const int startCol = groups * columnsPerLine + 1;
575  const int stopCol = min( (groups + 1) * columnsPerLine , L );
576  thestream << prefix << "Columns " << startCol << " to " << stopCol << "\n \n";
577  for (int row=0; row<L; row++){
578  for (int col=startCol-1; col<stopCol; col++){
579  if ((row==col) && (table==MutInfo)){
580  thestream << prefix << " 0 ";
581  for (int cnt=0; cnt<sPrecision; cnt++){ thestream << " "; }
582  } else {
583  int row2 = (Prob->gReorder()) ? Prob->gf1(row) : row;
584  int col2 = (Prob->gReorder()) ? Prob->gf1(col) : col;
585  if (table[row2 + L*col2] < 0.0){
586  thestream << prefix << table[row2 + L*col2];
587  } else {
588  thestream << prefix << " " << table[row2 + L*col2];
589  }
590  }
591  }
592  thestream << "\n";
593  }
594  thestream << "\n";
595 
596  }
597 
598  cout << thestream.str();
599 
600 }
601 
602 void CheMPS2::Correlations::Print(const int precision, const int columnsPerLine) const{
603 
604  cout << "--------------------------------------------------------" << endl;
605  cout << "Spin correlation function = 4 * ( < S_i^z S_j^z > - < S_i^z > * < S_j^z > ) \nHamiltonian index order is used!\n" << endl;
606  PrintTableNice( Cspin , precision, columnsPerLine );
607  cout << "--------------------------------------------------------" << endl;
608  cout << "Spin-flip correlation function = < S_i^+ S_j^- > + < S_i^- S_j^+ > \nHamiltonian index order is used!\n" << endl;
609  PrintTableNice( Cspinflip , precision, columnsPerLine);
610  cout << "--------------------------------------------------------" << endl;
611  cout << "Density correlation function = < n_i n_j > - < n_i > * < n_j > \nHamiltonian index order is used!\n" << endl;
612  PrintTableNice( Cdens , precision, columnsPerLine );
613  cout << "--------------------------------------------------------" << endl;
614  cout << "Singlet diradical correlation function = < d_i,up d_j,down > + < d_i,down d_j,up > - < d_i,up > * < d_j,down > - < d_i,down > * < d_j,up > \nHamiltonian index order is used!\n" << endl;
615  PrintTableNice( Cdirad , precision, columnsPerLine );
616  cout << "--------------------------------------------------------" << endl;
617  cout << "Two-orbital mutual information = 0.5 * ( s1(i) + s1(j) - s2(i,j) ) * ( 1 - delta(i,j) ) \nHamiltonian index order is used!\n" << endl;
618  PrintTableNice( MutInfo , precision, columnsPerLine );
619  cout << "--------------------------------------------------------" << endl;
620 
621 }
622 
623 
int gNmin(const int boundary) const
Get the min. possible particle number for a certain boundary.
Correlations(const SyBookkeeper *denBKIn, const Problem *ProbIn, TwoDM *the2DMin)
Constructor.
bool gReorder() const
Get whether the Hamiltonian orbitals are reordered for the DMRG calculation.
Definition: Problem.cpp:346
int gIndex() const
Get the location index.
Definition: TensorT.cpp:161
double getMutualInformation_HAM(const int row, const int col) const
Get a mutual information term, using the HAM indices.
double getMutualInformation_DMRG(const int row, const int col) const
Get a mutual information term, using the DMRG indices.
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.
void Print(const int precision=6, const int columnsPerLine=8) const
Print the correlation functions and two-orbital mutual information.
double getCdens_HAM(const int row, const int col) const
Get a Cdens term, using the HAM indices.
double getCdirad_HAM(const int row, const int col) const
Get a Cdirad term, using the HAM indices.
double SingleOrbitalEntropy_HAM(const int index) const
Get the single-orbital entropy for a certain site, using hte HAM indices.
double getTwoDMB_DMRG(const int cnt1, const int cnt2, const int cnt3, const int cnt4) const
Get a 2DM_B term, using the DMRG indices.
Definition: TwoDM.cpp:112
double SingleOrbitalEntropy_DMRG(const int index) const
Get the single-orbital entropy for a certain site, using the DMRG indices.
void mpi_broadcast()
Broadcast the diradical correlation function and the two-orbital mutual information.
double getCspin_HAM(const int row, const int col) const
Get a Cspin term, using the HAM indices.
int gTwoSmin(const int boundary, const int N) const
Get the minimum possible spin value for a certain boundary and particle number.
int gTwoSmax(const int boundary, const int N) const
Get the maximum possible spin value for a certain boundary and particle number.
double getCspinflip_HAM(const int row, const int col) const
Get a Cspinflip term, using the HAM indices.
double getCspin_DMRG(const int row, const int col) const
Get a Cspin term, using the DMRG indices.
int gTwoS() const
Get twice the targeted spin.
Definition: Problem.h:64
double * gStorage()
Get the pointer to the storage.
virtual ~Correlations()
Destructor.
double get1RDM_DMRG(const int cnt1, const int cnt2) const
Get a 1-RDM term, using the DMRG indices.
Definition: TwoDM.cpp:127
double getTwoDMA_DMRG(const int cnt1, const int cnt2, const int cnt3, const int cnt4) const
Get a 2DM_A term, using the DMRG indices.
Definition: TwoDM.cpp:97
double getCspinflip_DMRG(const int row, const int col) const
Get a Cspinflip term, using the DMRG indices.
double * gStorage()
Get the pointer to the storage.
Definition: TensorT.cpp:134
double MutualInformationDistance(const double power) const
Return Idistance(power) (see return for definition)
int gf1(const int HamOrb) const
Get the DMRG index corresponding to a Ham index.
Definition: Problem.cpp:347
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
double getCdens_DMRG(const int row, const int col) const
Get a Cdens term, using the DMRG indices.
int getNumberOfIrreps() const
Get the total number of irreps.
int gL() const
Get the number of orbitals.
int gMaxDimAtBound(const int boundary) const
Get the maximum virtual dimension at a certain boundary.
int gNmax(const int boundary) const
Get the max. possible particle number for a certain boundary.
int gIrrep(const int orbital) const
Get an orbital irrep.
double getCdirad_DMRG(const int row, const int col) const
Get a Cdirad term, using the DMRG indices.