CheMPS2
CASSCFpt2.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 <string>
23 #include <math.h>
24 #include <algorithm>
25 #include <sys/stat.h>
26 #include <assert.h>
27 
28 #include "CASSCF.h"
29 #include "DMRG.h"
30 #include "FCI.h"
31 #include "DMRGSCFrotations.h"
32 #include "Lapack.h"
33 #include "Cumulant.h"
34 #include "CASPT2.h"
35 #include "MPIchemps2.h"
36 #include "EdmistonRuedenberg.h"
37 
38 using std::string;
39 using std::ifstream;
40 using std::cout;
41 using std::endl;
42 using std::max;
43 
44 void CheMPS2::CASSCF::write_f4rdm_checkpoint( const string f4rdm_file, int * hamorb1, int * hamorb2, const int tot_dmrg_power6, double * contract ){
45 
46  #ifdef CHEMPS2_MPI_COMPILATION
47  const bool am_i_master = ( MPIchemps2::mpi_rank() == MPI_CHEMPS2_MASTER );
48  #else
49  const bool am_i_master = true;
50  #endif
51 
52  if ( am_i_master ){
53 
54  hid_t file_id = H5Fcreate( f4rdm_file.c_str(), H5F_ACC_TRUNC, H5P_DEFAULT, H5P_DEFAULT );
55  hid_t group_id = H5Gcreate( file_id, "/F4RDM", H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT );
56 
57  hsize_t dimarray1 = 1;
58  hid_t dataspace1_id = H5Screate_simple( 1, &dimarray1, NULL );
59  hid_t dataset1_id = H5Dcreate( group_id, "hamorb1", H5T_NATIVE_INT, dataspace1_id, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT );
60  H5Dwrite( dataset1_id, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, hamorb1 );
61  H5Dclose( dataset1_id );
62  H5Sclose( dataspace1_id );
63 
64  hsize_t dimarray2 = 1;
65  hid_t dataspace2_id = H5Screate_simple( 1, &dimarray2, NULL );
66  hid_t dataset2_id = H5Dcreate( group_id, "hamorb2", H5T_NATIVE_INT, dataspace2_id, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT );
67  H5Dwrite( dataset2_id, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, hamorb2 );
68  H5Dclose( dataset2_id );
69  H5Sclose( dataspace2_id );
70 
71  hsize_t dimarray3 = tot_dmrg_power6;
72  hid_t dataspace3_id = H5Screate_simple( 1, &dimarray3, NULL );
73  hid_t dataset3_id = H5Dcreate( group_id, "contract", H5T_NATIVE_DOUBLE, dataspace3_id, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT );
74  H5Dwrite( dataset3_id, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, H5P_DEFAULT, contract );
75  H5Dclose( dataset3_id );
76  H5Sclose( dataspace3_id );
77 
78  H5Gclose( group_id );
79  H5Fclose( file_id );
80 
81  cout << "Created F.4-RDM checkpoint file " << f4rdm_file << " at next orbitals ( " << hamorb1[ 0 ] << " , " << hamorb2[ 0 ] << " )." << endl;
82 
83  }
84 
85 }
86 
87 bool CheMPS2::CASSCF::read_f4rdm_checkpoint( const string f4rdm_file, int * hamorb1, int * hamorb2, const int tot_dmrg_power6, double * contract ){
88 
89  #ifdef CHEMPS2_MPI_COMPILATION
90  const bool am_i_master = ( MPIchemps2::mpi_rank() == MPI_CHEMPS2_MASTER );
91  #else
92  const bool am_i_master = true;
93  #endif
94 
95  // Check whether the file exists
96  int exists = 0;
97  if ( am_i_master ){
98  struct stat file_info;
99  const int file_stat = stat( f4rdm_file.c_str(), &file_info );
100  if ( file_stat == 0 ){ exists = 1; }
101  }
102  #ifdef CHEMPS2_MPI_COMPILATION
103  MPIchemps2::broadcast_array_int( &exists, 1, MPI_CHEMPS2_MASTER );
104  #endif
105  if ( exists == 0 ){
106  return false; // Not loaded
107  }
108 
109  if ( am_i_master ){
110 
111  hid_t file_id = H5Fopen( f4rdm_file.c_str(), H5F_ACC_RDONLY, H5P_DEFAULT );
112  hid_t group_id = H5Gopen( file_id, "/F4RDM", H5P_DEFAULT );
113 
114  hid_t dataset1_id = H5Dopen( group_id, "hamorb1", H5P_DEFAULT );
115  H5Dread( dataset1_id, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, hamorb1 );
116  H5Dclose( dataset1_id );
117 
118  hid_t dataset2_id = H5Dopen( group_id, "hamorb2", H5P_DEFAULT );
119  H5Dread( dataset2_id, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, hamorb2 );
120  H5Dclose( dataset2_id );
121 
122  hid_t dataset3_id = H5Dopen( group_id, "contract", H5P_DEFAULT );
123  H5Dread( dataset3_id, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, H5P_DEFAULT, contract );
124  H5Dclose( dataset3_id );
125 
126  H5Gclose( group_id );
127  H5Fclose( file_id );
128 
129  }
130  #ifdef CHEMPS2_MPI_COMPILATION
131  MPIchemps2::broadcast_array_int( hamorb1, 1, MPI_CHEMPS2_MASTER );
132  MPIchemps2::broadcast_array_int( hamorb2, 1, MPI_CHEMPS2_MASTER );
133  MPIchemps2::broadcast_array_double( contract, tot_dmrg_power6, MPI_CHEMPS2_MASTER );
134  #endif
135 
136  return true; // Loaded
137 
138 }
139 
140 void CheMPS2::CASSCF::fock_dot_4rdm( double * fockmx, CheMPS2::DMRG * dmrgsolver, CheMPS2::Hamiltonian * ham, int next_orb1, int next_orb2, double * work, double * result, const bool CHECKPOINT, const bool PSEUDOCANONICAL ){
141 
142  const int LAS = ham->getL();
143  int size = LAS * LAS * LAS * LAS * LAS * LAS;
144  int inc1 = 1;
145 
146  for ( int diag = 0; diag < LAS; diag++ ){
147  if (( next_orb1 == diag ) && ( next_orb2 == diag )){
148  double prefactor = 0.5 * fockmx[ diag + LAS * diag ];
149  if ( fabs( prefactor ) > 0.0 ){
150  dmrgsolver->Symm4RDM( work, diag, diag, false );
151  daxpy_( &size, &prefactor, work, &inc1, result, &inc1 );
152  }
153  if ( diag == LAS - 1 ){
154  next_orb1 = 0;
155  next_orb2 = 1;
156  } else {
157  next_orb1 = diag + 1;
158  next_orb2 = diag + 1;
159  }
160  if ( CHECKPOINT ){ write_f4rdm_checkpoint( CheMPS2::DMRGSCF_f4rdm_name, &next_orb1, &next_orb2, size, result ); }
161  }
162  }
163 
164  if ( PSEUDOCANONICAL == false ){
165  for ( int orb1 = 0; orb1 < LAS; orb1++ ){
166  for ( int orb2 = orb1 + 1; orb2 < LAS; orb2++ ){
167  if (( next_orb1 == orb1 ) && ( next_orb2 == orb2 )){
168  double prefactor = 0.5 * ( fockmx[ orb1 + LAS * orb2 ] + fockmx[ orb2 + LAS * orb1 ] );
169  if (( ham->getOrbitalIrrep( orb1 ) == ham->getOrbitalIrrep( orb2 ) ) && ( fabs( prefactor ) > 0.0 )){
170  dmrgsolver->Symm4RDM( work, orb1, orb2, false );
171  daxpy_( &size, &prefactor, work, &inc1, result, &inc1 );
172  }
173  if ( orb2 == LAS - 1 ){
174  next_orb1 = next_orb1 + 1;
175  next_orb2 = next_orb1 + 1;
176  } else {
177  next_orb2 = next_orb2 + 1;
178  }
179  if (( ham->getOrbitalIrrep( orb1 ) == ham->getOrbitalIrrep( orb2 ) ) && ( CHECKPOINT )){
180  write_f4rdm_checkpoint( CheMPS2::DMRGSCF_f4rdm_name, &next_orb1, &next_orb2, size, result );
181  }
182  }
183  }
184  }
185  }
186 
187 }
188 
189 double CheMPS2::CASSCF::caspt2( const int Nelectrons, const int TwoS, const int Irrep, ConvergenceScheme * OptScheme, const int rootNum, DMRGSCFoptions * scf_options, const double IPEA, const double IMAG, const bool PSEUDOCANONICAL, const bool CHECKPOINT, const bool CUMULANT ){
190 
191  #ifdef CHEMPS2_MPI_COMPILATION
192  const bool am_i_master = ( MPIchemps2::mpi_rank() == MPI_CHEMPS2_MASTER );
193  #else
194  const bool am_i_master = true;
195  #endif
196 
197  const int num_elec = Nelectrons - 2 * iHandler->getNOCCsum();
198  assert( num_elec >= 0 );
199 
200  if ( CASPT2::vector_length( iHandler ) == 0 ){
201  if ( am_i_master ){
202  cout << "CheMPS2::CASSCF::caspt2 : There are no CASPT2 excitations between the CORE, ACTIVE, and VIRTUAL orbital spaces." << endl;
203  }
204  return 0.0;
205  }
206 
207  //Determine the maximum NORB(irrep) and the max_block_size for the ERI orbital rotation
208  const int maxlinsize = iHandler->getNORBmax();
209  const long long fullsize = ((long long) maxlinsize ) * ((long long) maxlinsize ) * ((long long) maxlinsize ) * ((long long) maxlinsize );
210  const string tmp_filename = tmp_folder + "/" + CheMPS2::DMRGSCF_eri_storage_name;
211  const int dmrgsize_power4 = nOrbDMRG * nOrbDMRG * nOrbDMRG * nOrbDMRG;
212  //For (ERI rotation, update unitary, block diagonalize, orbital localization)
213  DMRGSCFintegrals * theRotatedTEI = new DMRGSCFintegrals( iHandler );
214  const int temp_work_size = (( fullsize > CheMPS2::DMRGSCF_max_mem_eri_tfo ) ? CheMPS2::DMRGSCF_max_mem_eri_tfo : fullsize );
215  const int work_mem_size = max( max( temp_work_size , maxlinsize * maxlinsize * 4 ) , dmrgsize_power4 );
216  const int tot_dmrg_power6 = dmrgsize_power4 * nOrbDMRG * nOrbDMRG;
217  double * mem1 = new double[ work_mem_size ];
218  double * mem2 = new double[ ( PSEUDOCANONICAL ) ? work_mem_size : max( work_mem_size, tot_dmrg_power6 ) ];
219 
220  // If you did not run CheMPS2::CASSCF::solve, you NEED to load the unitary from disk
221  if ( successful_solve == false ){
222  assert( scf_options->getStoreUnitary() );
223  if ( am_i_master ){
224  struct stat file_info;
225  const int file_stat = stat( (scf_options->getUnitaryStorageName()).c_str(), &file_info );
226  assert( file_stat == 0 );
227  unitary->loadU( scf_options->getUnitaryStorageName() );
228  }
229  #ifdef CHEMPS2_MPI_COMPILATION
230  unitary->broadcast( MPI_CHEMPS2_MASTER );
231  #endif
232  } else {
233  // Rotate to pseudocanonical orbitals
234  if ( PSEUDOCANONICAL ){ // successful_solve needs to be true, because the 1-RDM is needed for buildQmatACT();
235  buildTmatrix();
236  buildQmatOCC();
237  buildQmatACT();
238  construct_fock( theFmatrix, theTmatrix, theQmatOCC, theQmatACT, iHandler );
239  block_diagonalize( 'O', theFmatrix, unitary, mem1, mem2, iHandler, false, NULL, NULL, NULL );
240  block_diagonalize( 'A', theFmatrix, unitary, mem1, mem2, iHandler, false, NULL, NULL, NULL );
241  block_diagonalize( 'V', theFmatrix, unitary, mem1, mem2, iHandler, false, NULL, NULL, NULL );
242  if (( am_i_master ) && ( scf_options->getStoreUnitary() )){ unitary->saveU( scf_options->getUnitaryStorageName() ); }
243  }
244  }
245 
246  // Fill active space Hamiltonian
247  Hamiltonian * HamAS = new Hamiltonian( nOrbDMRG, SymmInfo.getGroupNumber(), iHandler->getIrrepOfEachDMRGorbital() );
248  Problem * Prob = new Problem( HamAS, TwoS, num_elec, Irrep );
249  Prob->SetupReorderD2h(); // Doesn't matter if the group isn't D2h, Prob checks it.
250  buildTmatrix();
251  buildQmatOCC();
252  fillConstAndTmatDMRG( HamAS );
253  if ( am_i_master ){
254  DMRGSCFrotations::rotate( VMAT_ORIG, HamAS->getVmat(), NULL, 'A', 'A', 'A', 'A', iHandler, unitary, mem1, mem2, work_mem_size, tmp_filename );
255  }
256  #ifdef CHEMPS2_MPI_COMPILATION
257  HamAS->getVmat()->broadcast( MPI_CHEMPS2_MASTER );
258  #endif
259 
260  // Reorder the orbitals based on the Fiedler vector of the exchange matrix
261  if ( scf_options->getWhichActiveSpace() == 3 ){
262  int * dmrg2ham = new int[ nOrbDMRG ];
263  if ( am_i_master ){
264  EdmistonRuedenberg * theLocalizer = new EdmistonRuedenberg( HamAS->getVmat(), iHandler->getGroupNumber() );
265  theLocalizer->FiedlerGlobal( dmrg2ham );
266  delete theLocalizer;
267  }
268  #ifdef CHEMPS2_MPI_COMPILATION
269  MPIchemps2::broadcast_array_int( dmrg2ham, nOrbDMRG, MPI_CHEMPS2_MASTER );
270  #endif
271  Prob->setup_reorder_custom( dmrg2ham );
272  delete [] dmrg2ham;
273  }
274 
275  double E_CASSCF = 0.0;
276  double * three_dm = new double[ tot_dmrg_power6 ];
277  double * contract = new double[ tot_dmrg_power6 ];
278  for ( int cnt = 0; cnt < tot_dmrg_power6; cnt++ ){ contract[ cnt ] = 0.0; }
279 
280  int next_hamorb1 = 0;
281  int next_hamorb2 = 0;
282  const bool make_checkpt = (( CUMULANT == false ) && ( CHECKPOINT ));
283  bool checkpt_loaded = false;
284  if ( make_checkpt ){
285  assert(( OptScheme != NULL ) || ( rootNum > 1 ));
286  checkpt_loaded = read_f4rdm_checkpoint( CheMPS2::DMRGSCF_f4rdm_name, &next_hamorb1, &next_hamorb2, tot_dmrg_power6, contract );
287  }
288 
289  // Solve the active space problem
290  if (( OptScheme == NULL ) && ( rootNum == 1 )){ // Do FCI
291 
292  if ( am_i_master ){
293  const int nalpha = ( num_elec + TwoS ) / 2;
294  const int nbeta = ( num_elec - TwoS ) / 2;
295  const double workmem = 1000.0; // 1GB
296  const int verbose = 2;
297  CheMPS2::FCI * theFCI = new CheMPS2::FCI( HamAS, nalpha, nbeta, Irrep, workmem, verbose );
298  double * inoutput = new double[ theFCI->getVecLength(0) ];
299  theFCI->ClearVector( theFCI->getVecLength(0), inoutput );
300  inoutput[ theFCI->LowestEnergyDeterminant() ] = 1.0;
301  E_CASSCF = theFCI->GSDavidson( inoutput );
302  theFCI->Fill2RDM( inoutput, DMRG2DM ); // 2-RDM
303  theFCI->Fill3RDM( inoutput, three_dm ); // 3-RDM
304  setDMRG1DM( num_elec, nOrbDMRG, DMRG1DM, DMRG2DM ); // 1-RDM
305  buildQmatACT();
306  construct_fock( theFmatrix, theTmatrix, theQmatOCC, theQmatACT, iHandler );
307  copy_active( theFmatrix, mem2, iHandler ); // Fock
308  theFCI->Fock4RDM( inoutput, three_dm, mem2, contract ); // trace( Fock * 4-RDM )
309  delete theFCI;
310  delete [] inoutput;
311  }
312  #ifdef CHEMPS2_MPI_COMPILATION
313  MPIchemps2::broadcast_array_double( &E_CASSCF, 1, MPI_CHEMPS2_MASTER );
314  MPIchemps2::broadcast_array_double( DMRG2DM, dmrgsize_power4, MPI_CHEMPS2_MASTER );
315  MPIchemps2::broadcast_array_double( three_dm, tot_dmrg_power6, MPI_CHEMPS2_MASTER );
316  MPIchemps2::broadcast_array_double( contract, tot_dmrg_power6, MPI_CHEMPS2_MASTER );
317  setDMRG1DM( num_elec, nOrbDMRG, DMRG1DM, DMRG2DM );
318  #endif
319 
320  } else { // Do the DMRG sweeps
321 
322  assert( OptScheme != NULL );
323  for ( int cnt = 0; cnt < dmrgsize_power4; cnt++ ){ DMRG2DM[ cnt ] = 0.0; } // Clear the 2-RDM
324  CheMPS2::DMRG * theDMRG = new DMRG( Prob, OptScheme, make_checkpt, tmp_folder );
325  for ( int state = 0; state < rootNum; state++ ){
326  if ( state > 0 ){ theDMRG->newExcitation( fabs( E_CASSCF ) ); }
327  if ( checkpt_loaded == false ){ E_CASSCF = theDMRG->Solve(); }
328  if (( state == 0 ) && ( rootNum > 1 )){ theDMRG->activateExcitations( rootNum - 1 ); }
329  }
330  theDMRG->calc_rdms_and_correlations( true );
331  copy2DMover( theDMRG->get2DM(), nOrbDMRG, DMRG2DM ); // 2-RDM
332  setDMRG1DM( num_elec, nOrbDMRG, DMRG1DM, DMRG2DM ); // 1-RDM
333  buildQmatACT();
334  construct_fock( theFmatrix, theTmatrix, theQmatOCC, theQmatACT, iHandler );
335  copy_active( theFmatrix, mem2, iHandler ); // Fock
336  if ( CUMULANT ){
337  CheMPS2::Cumulant::gamma4_fock_contract_ham( Prob, theDMRG->get3DM(), theDMRG->get2DM(), mem2, contract );
338  } else {
339  fock_dot_4rdm( mem2, theDMRG, HamAS, next_hamorb1, next_hamorb2, three_dm, contract, make_checkpt, PSEUDOCANONICAL );
340  }
341  theDMRG->get3DM()->fill_ham_index( 1.0, false, three_dm, 0, nOrbDMRG );
342  if (( CheMPS2::DMRG_storeMpsOnDisk ) && ( make_checkpt == false )){ theDMRG->deleteStoredMPS(); }
343  if ( CheMPS2::DMRG_storeRenormOptrOnDisk ){ theDMRG->deleteStoredOperators(); }
344  delete theDMRG;
345 
346  }
347 
348  delete Prob;
349  delete HamAS;
350 
351  if ( PSEUDOCANONICAL == false ){
352  if ( am_i_master ){ cout << "CASPT2 : Deviation from pseudocanonical = " << deviation_from_blockdiag( theFmatrix, iHandler ) << endl; }
353  block_diagonalize( 'O', theFmatrix, unitary, mem1, mem2, iHandler, false, NULL, NULL, NULL );
354  block_diagonalize( 'A', theFmatrix, unitary, mem1, mem2, iHandler, false, DMRG2DM, three_dm, contract ); // 2-RDM, 3-RDM, and trace( Fock * cu(4)-4-RDM )
355  block_diagonalize( 'V', theFmatrix, unitary, mem1, mem2, iHandler, false, NULL, NULL, NULL );
356  setDMRG1DM( num_elec, nOrbDMRG, DMRG1DM, DMRG2DM ); // 1-RDM
357  buildTmatrix();
358  buildQmatOCC();
359  buildQmatACT();
360  construct_fock( theFmatrix, theTmatrix, theQmatOCC, theQmatACT, iHandler ); // Fock
361  }
362 
363  // Calculate the matrix elements needed to calculate the CASPT2 V-vector
364  if ( am_i_master ){
365  DMRGSCFrotations::rotate( VMAT_ORIG, NULL, theRotatedTEI, 'C', 'C', 'F', 'F', iHandler, unitary, mem1, mem2, work_mem_size, tmp_filename );
366  DMRGSCFrotations::rotate( VMAT_ORIG, NULL, theRotatedTEI, 'C', 'V', 'C', 'V', iHandler, unitary, mem1, mem2, work_mem_size, tmp_filename );
367  delete_file( tmp_filename );
368  }
369 
370  delete [] mem1;
371  delete [] mem2;
372 
373  double E_CASPT2 = 0.0;
374  if ( am_i_master ){
375  cout << "CASPT2 : Deviation from pseudocanonical = " << deviation_from_blockdiag( theFmatrix, iHandler ) << endl;
376  CheMPS2::CASPT2 * myCASPT2 = new CheMPS2::CASPT2( iHandler, theRotatedTEI, theTmatrix, theFmatrix, DMRG1DM, DMRG2DM, three_dm, contract, IPEA );
377  delete theRotatedTEI;
378  delete [] three_dm;
379  delete [] contract;
380  E_CASPT2 = myCASPT2->solve( IMAG );
381  delete myCASPT2;
382  } else {
383  delete theRotatedTEI;
384  delete [] three_dm;
385  delete [] contract;
386  }
387  #ifdef CHEMPS2_MPI_COMPILATION
388  MPIchemps2::broadcast_array_double( &E_CASPT2, 1, MPI_CHEMPS2_MASTER );
389  #endif
390 
391  return E_CASPT2;
392 
393 }
394 
395 void CheMPS2::CASSCF::construct_fock( DMRGSCFmatrix * Fock, const DMRGSCFmatrix * Tmat, const DMRGSCFmatrix * Qocc, const DMRGSCFmatrix * Qact, const DMRGSCFindices * idx ){
396 
397  const int n_irreps = idx->getNirreps();
398  for ( int irrep = 0; irrep < n_irreps; irrep++ ){
399  const int NORB = idx->getNORB( irrep );
400  for (int row = 0; row < NORB; row++){
401  for (int col = 0; col < NORB; col++){
402  Fock->set( irrep, row, col, Tmat->get( irrep, row, col )
403  + Qocc->get( irrep, row, col )
404  + Qact->get( irrep, row, col ) );
405  }
406  }
407  }
408 
409 }
410 
412 
413  double rms_deviation = 0.0;
414  const int n_irreps = idx->getNirreps();
415  for ( int irrep = 0; irrep < n_irreps; irrep++ ){
416  const int NOCC = idx->getNOCC( irrep );
417  for ( int row = 0; row < NOCC; row++ ){
418  for ( int col = row + 1; col < NOCC; col++ ){
419  rms_deviation += matrix->get( irrep, row, col ) * matrix->get( irrep, row, col );
420  rms_deviation += matrix->get( irrep, col, row ) * matrix->get( irrep, col, row );
421  }
422  }
423  const int NACT = idx->getNDMRG( irrep );
424  for ( int row = 0; row < NACT; row++ ){
425  for ( int col = row + 1; col < NACT; col++ ){
426  rms_deviation += matrix->get( irrep, NOCC + row, NOCC + col ) * matrix->get( irrep, NOCC + row, NOCC + col );
427  rms_deviation += matrix->get( irrep, NOCC + col, NOCC + row ) * matrix->get( irrep, NOCC + col, NOCC + row );
428  }
429  }
430  const int NVIR = idx->getNVIRT( irrep );
431  const int N_OA = NOCC + NACT;
432  for ( int row = 0; row < NVIR; row++ ){
433  for ( int col = row + 1; col < NVIR; col++ ){
434  rms_deviation += matrix->get( irrep, N_OA + row, N_OA + col ) * matrix->get( irrep, N_OA + row, N_OA + col );
435  rms_deviation += matrix->get( irrep, N_OA + col, N_OA + row ) * matrix->get( irrep, N_OA + col, N_OA + row );
436  }
437  }
438  }
439  rms_deviation = sqrt( rms_deviation );
440  return rms_deviation;
441 
442 }
443 
444 
static void write_f4rdm_checkpoint(const string f4rdm_file, int *hamorb1, int *hamorb2, const int tot_dmrg_power6, double *contract)
Write the checkpoint file for the contraction of the generalized Fock operator with the 4-RDM to disk...
Definition: CASSCFpt2.cpp:44
int getNORB(const int irrep) const
Get the number of orbitals for an irrep.
int getNOCCsum() const
Get the total number of occupied orbitals.
static void rotate(const FourIndex *ORIG_VMAT, FourIndex *NEW_VMAT, DMRGSCFintegrals *ROT_TEI, const char space1, const char space2, const char space3, const char space4, DMRGSCFindices *idx, DMRGSCFunitary *umat, double *mem1, double *mem2, const int mem_size, const string filename)
Fill the rotated two-body matrix elements for the space. If the blocks become too large...
void Symm4RDM(double *output, const int ham_orb1, const int ham_orb2, const bool last_case)
Obtain the symmetrized 4-RDM terms 0.5 * ( Gamma4_ijkl,pqrt + Gamma4_ijkt,pqrl ) with l and t fixed...
Definition: DMRGfock.cpp:41
double GSDavidson(double *inoutput=NULL, const int DVDSN_NUM_VEC=CheMPS2::DAVIDSON_NUM_VEC) const
Calculates the FCI ground state with Davidson&#39;s algorithm.
Definition: FCI.cpp:1920
double solve(const double imag_shift, const bool CONJUGATE_GRADIENT=false) const
Solve for the CASPT2 energy (note that the IPEA shift has been set in the constructor) ...
Definition: CASPT2.cpp:206
static void copy_active(double *origin, DMRGSCFmatrix *result, const DMRGSCFindices *idx, const bool one_rdm)
Copy a one-orbital quantity from array format to DMRGSCFmatrix format.
Definition: CASSCF.cpp:375
void deleteStoredMPS()
Call "rm " + CheMPS2::DMRG_MPS_storage_prefix + "*.h5".
Definition: DMRGmpsio.cpp:164
int getGroupNumber() const
Get the group number.
static void construct_fock(DMRGSCFmatrix *Fock, const DMRGSCFmatrix *Tmat, const DMRGSCFmatrix *Qocc, const DMRGSCFmatrix *Qact, const DMRGSCFindices *idx)
Construct the Fock matrix.
Definition: CASSCFpt2.cpp:395
void loadU(const string filename=DMRGSCF_unitary_storage_name)
Load the unitary from disk.
int * getIrrepOfEachDMRGorbital()
Get an array with the irreps of each DMRG orbital.
int getOrbitalIrrep(const int nOrb) const
Get an orbital irrep number.
Definition: Hamiltonian.cpp:97
static int mpi_rank()
Get the rank of this MPI process.
Definition: MPIchemps2.h:131
static long long vector_length(const DMRGSCFindices *idx)
Return the vector length for the CASPT2 first order wavefunction (before diagonalization of the overl...
Definition: CASPT2.cpp:793
double getVmat(const int index1, const int index2, const int index3, const int index4) const
Get a Vmat element.
TwoDM * get2DM()
Get the pointer to the 2-RDM.
Definition: DMRG.h:113
int getGroupNumber() const
Get the group number.
Definition: Irreps.cpp:66
void set(const int irrep, const int p, const int q, const double val)
Set an element.
int getNORBmax() const
Get the maximum NORB.
unsigned int LowestEnergyDeterminant() const
Return the global counter of the Slater determinant with the lowest energy.
Definition: FCI.cpp:1700
void saveU(const string filename=DMRGSCF_unitary_storage_name) const
Save the unitary to disk.
void fill_ham_index(const double alpha, const bool add, double *storage, const int last_orb_start, const int last_orb_num)
Perform storage[ :, :, :, :, :, : ] { = or += } alpha * 3-RDM[ :, :, :, :, :, last_orb_start: last_or...
Definition: ThreeDM.cpp:155
double Solve()
Solver.
Definition: DMRG.cpp:199
static void copy2DMover(TwoDM *theDMRG2DM, const int LAS, double *two_dm)
Copy over the DMRG 2-RDM.
Definition: CASSCF.cpp:119
int getWhichActiveSpace() const
Get which active space should be considered in the DMRG routine.
double caspt2(const int Nelectrons, const int TwoS, const int Irrep, ConvergenceScheme *OptScheme, const int rootNum, DMRGSCFoptions *scf_options, const double IPEA, const double IMAG, const bool PSEUDOCANONICAL, const bool CHECKPOINT=false, const bool CUMULANT=false)
Calculate the caspt2 correction energy for a converged casscf wavefunction.
Definition: CASSCFpt2.cpp:189
static void block_diagonalize(const char space, const DMRGSCFmatrix *Mat, DMRGSCFunitary *Umat, double *work1, double *work2, const DMRGSCFindices *idx, const bool invert, double *two_dm, double *three_dm, double *contract)
Block-diagonalize Mat.
Definition: CASSCF.cpp:433
void FiedlerGlobal(int *dmrg2ham) const
Permute the orbitals so that the bandwidth of the exchange matrix is (approximately) minimized...
void newExcitation(const double EshiftIn)
Push back current calculation and set everything up to calculate a (new) excitation.
Definition: DMRG.cpp:395
static void fock_dot_4rdm(double *fockmx, CheMPS2::DMRG *dmrgsolver, CheMPS2::Hamiltonian *ham, int next_orb1, int next_orb2, double *work, double *result, const bool CHECKPOINT, const bool PSEUDOCANONICAL)
Build the contraction of the fock matrix with the 4-RDM.
Definition: CASSCFpt2.cpp:140
double get(const int irrep, const int p, const int q) const
Get an element.
bool getStoreUnitary() const
Get whether the Orbital Rotation checkpoint should be stored to disk.
unsigned int getVecLength(const int irrep_center) const
Getter for the number of variables in the vector " E_ij | FCI vector > " ; where irrep_center = I_i x...
Definition: FCI.h:65
static void setDMRG1DM(const int num_elec, const int LAS, double *one_dm, double *two_dm)
Construct the 1-RDM from the 2-RDM.
Definition: CASSCF.cpp:134
void Fill3RDM(double *vector, double *ThreeRDM) const
Construct the (spin-summed) 3-RDM of a FCI vector: Gamma^3(i,j,k,l,m,n) = sum_sigma,tau,s < a^+_{i,sigma} a^+_{j,tau} a^+_{k,s} a_{n,s} a_{m,tau} a_{l,sigma} > = ThreeRDM[ i + L * ( j + L * ( k + L * ( l + L * ( m + L * n ) ) ) ) ].
Definition: FCI.cpp:1168
static void ClearVector(const unsigned int vecLength, double *vec)
Set the entries of a vector to zero.
Definition: FCI.cpp:1908
void activateExcitations(const int maxExcIn)
Activate the necessary storage and machinery to handle excitations.
Definition: DMRG.cpp:384
static void gamma4_fock_contract_ham(const Problem *prob, const ThreeDM *the3DM, const TwoDM *the2DM, double *fock, double *result)
Contract the CASPT2 Fock operator with the cumulant approximation of in time, using HAM indices...
Definition: Cumulant.cpp:84
int getNDMRG(const int irrep) const
Get the number of active orbitals for an irrep.
ThreeDM * get3DM()
Get the pointer to the 3-RDM.
Definition: DMRG.h:117
int getNOCC(const int irrep) const
Get the number of occupied orbitals for an irrep.
static bool read_f4rdm_checkpoint(const string f4rdm_file, int *hamorb1, int *hamorb2, const int tot_dmrg_power6, double *contract)
Read the checkpoint file for the contraction of the generalized Fock operator with the 4-RDM from dis...
Definition: CASSCFpt2.cpp:87
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...
void deleteStoredOperators()
Call "rm " + tempfolder + "/" + CheMPS2::DMRG_OPERATOR_storage_prefix + string(thePID) + "*...
static double deviation_from_blockdiag(DMRGSCFmatrix *matrix, const DMRGSCFindices *idx)
Return the RMS deviation from block-diagonal.
Definition: CASSCFpt2.cpp:411
double Fill2RDM(double *vector, double *TwoRDM) const
Construct the (spin-summed) 2-RDM of a FCI vector: Gamma^2(i,j,k,l) = sum_sigma,tau < a^+_i...
Definition: FCI.cpp:805
int getNVIRT(const int irrep) const
Get the number of virtual orbitals for an irrep.
void Fock4RDM(double *vector, double *ThreeRDM, double *Fock, double *output) const
Construct the (spin-summed) contraction of the 4-RDM with the Fock operator: output(i,j,k,p,q,r) = sum_{l,t} Fock(l,t) * Gamma^4(i,j,k,l,p,q,r,t)
Definition: FCI.cpp:1160
int getL() const
Get the number of orbitals.
Definition: Hamiltonian.cpp:93
string getUnitaryStorageName() const
Get the Orbital Rotation checkpoint filename.
int getNirreps() const
Get the number of irreps.