CheMPS2
Hamiltonian.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 <stdio.h>
21 #include <stdlib.h>
22 #include <assert.h>
23 #include <iostream>
24 #include <string>
25 #include <fstream>
26 #include <sys/stat.h>
27 
28 #include "Irreps.h"
29 #include "TwoIndex.h"
30 #include "FourIndex.h"
31 #include "Hamiltonian.h"
32 #include "MyHDF5.h"
33 
34 using std::cout;
35 using std::endl;
36 using std::string;
37 using std::ifstream;
38 
39 CheMPS2::Hamiltonian::Hamiltonian(const int Norbitals, const int nGroup, const int * OrbIrreps){
40 
41  L = Norbitals;
42  assert( nGroup>=0 );
43  assert( nGroup<=7 );
44  SymmInfo.setGroup(nGroup);
45 
46  orb2irrep = new int[L];
47  orb2indexSy = new int[L];
48  int nIrreps = SymmInfo.getNumberOfIrreps();
49  irrep2num_orb = new int[nIrreps];
50  for (int cnt=0; cnt<nIrreps; cnt++) irrep2num_orb[cnt] = 0;
51  for (int cnt=0; cnt<L; cnt++){
52  assert( OrbIrreps[cnt]>=0 );
53  assert( OrbIrreps[cnt]<nIrreps );
54  orb2irrep[cnt] = OrbIrreps[cnt];
55  orb2indexSy[cnt] = irrep2num_orb[orb2irrep[cnt]];
56  irrep2num_orb[orb2irrep[cnt]]++;
57  }
58 
59  Econst = 0.0;
60  Tmat = new TwoIndex(SymmInfo.getGroupNumber(),irrep2num_orb);
61  Vmat = new FourIndex(SymmInfo.getGroupNumber(),irrep2num_orb);
62 
63 }
64 
65 CheMPS2::Hamiltonian::Hamiltonian( const string filename, const int psi4groupnumber ){
66 
67  SymmInfo.setGroup( psi4groupnumber );
68  CreateAndFillFromFCIDUMP( filename );
69 
70 }
71 
72 CheMPS2::Hamiltonian::Hamiltonian( const bool fileh5, const string main_file, const string file_tmat, const string file_vmat ){
73 
74  if ( fileh5 ){
75  CreateAndFillFromH5( main_file, file_tmat, file_vmat );
76  } else {
77  cout << "CheMPS2::Hamiltonian::Hamiltonian( false, const string , const string , const string ) was deprecated." << endl;
78  assert( fileh5 == true );
79  }
80 
81 }
82 
84 
85  delete [] orb2irrep;
86  delete [] orb2indexSy;
87  delete [] irrep2num_orb;
88  delete Tmat;
89  delete Vmat;
90 
91 }
92 
93 int CheMPS2::Hamiltonian::getL() const{ return L; }
94 
95 int CheMPS2::Hamiltonian::getNGroup() const{ return SymmInfo.getGroupNumber(); }
96 
97 int CheMPS2::Hamiltonian::getOrbitalIrrep(const int nOrb) const{ return orb2irrep[nOrb]; }
98 
99 void CheMPS2::Hamiltonian::setEconst(const double val){ Econst = val; }
100 
101 double CheMPS2::Hamiltonian::getEconst() const{ return Econst; }
102 
103 void CheMPS2::Hamiltonian::setTmat(const int index1, const int index2, const double val){
104 
105  assert( orb2irrep[index1]==orb2irrep[index2] );
106  Tmat->set(orb2irrep[index1], orb2indexSy[index1], orb2indexSy[index2], val);
107 
108 }
109 
110 double CheMPS2::Hamiltonian::getTmat(const int index1, const int index2) const{
111 
112  if (orb2irrep[index1]==orb2irrep[index2]){
113  return Tmat->get(orb2irrep[index1], orb2indexSy[index1], orb2indexSy[index2]);
114  }
115 
116  return 0.0;
117 
118 }
119 
121 
122 void CheMPS2::Hamiltonian::setVmat(const int index1, const int index2, const int index3, const int index4, const double val){
123 
124  assert( Irreps::directProd(orb2irrep[index1],orb2irrep[index2]) == Irreps::directProd(orb2irrep[index3],orb2irrep[index4]) );
125  Vmat->set(orb2irrep[index1], orb2irrep[index2], orb2irrep[index3], orb2irrep[index4], orb2indexSy[index1], orb2indexSy[index2], orb2indexSy[index3], orb2indexSy[index4], val);
126 
127 }
128 
129 void CheMPS2::Hamiltonian::addToVmat(const int index1, const int index2, const int index3, const int index4, const double val){
130 
131  assert( Irreps::directProd(orb2irrep[index1],orb2irrep[index2]) == Irreps::directProd(orb2irrep[index3],orb2irrep[index4]) );
132  Vmat->add(orb2irrep[index1], orb2irrep[index2], orb2irrep[index3], orb2irrep[index4], orb2indexSy[index1], orb2indexSy[index2], orb2indexSy[index3], orb2indexSy[index4], val);
133 
134 }
135 
136 double CheMPS2::Hamiltonian::getVmat(const int index1, const int index2, const int index3, const int index4) const{
137 
138  if ( Irreps::directProd(orb2irrep[index1],orb2irrep[index2]) == Irreps::directProd(orb2irrep[index3],orb2irrep[index4]) ){
139  return Vmat->get(orb2irrep[index1], orb2irrep[index2], orb2irrep[index3], orb2irrep[index4], orb2indexSy[index1], orb2indexSy[index2], orb2indexSy[index3], orb2indexSy[index4]);
140  }
141 
142  return 0.0;
143 
144 }
145 
147 
148 void CheMPS2::Hamiltonian::save(const string file_parent, const string file_tmat, const string file_vmat) const{
149 
150  Tmat->save(file_tmat);
151  Vmat->save(file_vmat);
152 
153  //The hdf5 file
154  hid_t file_id = H5Fcreate(file_parent.c_str(), H5F_ACC_TRUNC, H5P_DEFAULT, H5P_DEFAULT);
155 
156  //The data
157  hid_t group_id = H5Gcreate(file_id, "/Data", H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
158 
159  //The chain length
160  hsize_t dimarray = 1;
161  hid_t dataspace_id = H5Screate_simple(1, &dimarray, NULL);
162  hid_t dataset_id = H5Dcreate(group_id, "L", H5T_STD_I32LE, dataspace_id, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
163  H5Dwrite(dataset_id, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &L);
164 
165  //The group number
166  hsize_t dimarray2 = 1;
167  hid_t dataspace_id2 = H5Screate_simple(1, &dimarray2, NULL);
168  hid_t dataset_id2 = H5Dcreate(group_id, "nGroup", H5T_STD_I32LE, dataspace_id2, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
169  int nGroup = SymmInfo.getGroupNumber();
170  H5Dwrite(dataset_id2, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &nGroup);
171 
172  //orb2irrep
173  hsize_t dimarray3 = L;
174  hid_t dataspace_id3 = H5Screate_simple(1, &dimarray3, NULL);
175  hid_t dataset_id3 = H5Dcreate(group_id, "orb2irrep", H5T_STD_I32LE, dataspace_id3, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
176  H5Dwrite(dataset_id3, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, orb2irrep);
177 
178  //Econst
179  hsize_t dimarray4 = 1;
180  hid_t dataspace_id4 = H5Screate_simple(1, &dimarray4, NULL);
181  hid_t dataset_id4 = H5Dcreate(group_id, "Econst", H5T_IEEE_F64LE, dataspace_id4, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
182  H5Dwrite(dataset_id4, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, H5P_DEFAULT, &Econst);
183 
184  H5Dclose(dataset_id);
185  H5Sclose(dataspace_id);
186  H5Dclose(dataset_id2);
187  H5Sclose(dataspace_id2);
188  H5Dclose(dataset_id3);
189  H5Sclose(dataspace_id3);
190  H5Dclose(dataset_id4);
191  H5Sclose(dataspace_id4);
192 
193  H5Gclose(group_id);
194 
195  H5Fclose(file_id);
196 
197 }
198 
199 void CheMPS2::Hamiltonian::read(const string file_parent, const string file_tmat, const string file_vmat){
200 
201  Tmat->read(file_tmat);
202  Vmat->read(file_vmat);
203 
204  //The hdf5 file
205  hid_t file_id = H5Fopen(file_parent.c_str(), H5F_ACC_RDONLY, H5P_DEFAULT);
206 
207  //The data
208  hid_t group_id = H5Gopen(file_id, "/Data",H5P_DEFAULT);
209 
210  //The chain length
211  hid_t dataset_id1 = H5Dopen(group_id, "L", H5P_DEFAULT);
212  int Lagain;
213  H5Dread(dataset_id1, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &Lagain);
214  assert( L==Lagain );
215 
216  //The group number
217  hid_t dataset_id2 = H5Dopen(group_id, "nGroup", H5P_DEFAULT);
218  int nGroup;
219  H5Dread(dataset_id2, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &nGroup);
220  assert( nGroup==SymmInfo.getGroupNumber() );
221 
222  //orb2irrep
223  hid_t dataset_id3 = H5Dopen(group_id, "orb2irrep", H5P_DEFAULT);
224  int * orb2irrepAgain = new int[L];
225  H5Dread(dataset_id3, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, orb2irrepAgain);
226  for (int cnt=0; cnt<L; cnt++){
227  assert( orb2irrep[cnt]==orb2irrepAgain[cnt] );
228  }
229  delete [] orb2irrepAgain;
230 
231  //Econst
232  hid_t dataset_id4 = H5Dopen(group_id, "Econst", H5P_DEFAULT);
233  H5Dread(dataset_id4, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, H5P_DEFAULT, &Econst);
234 
235  H5Dclose(dataset_id1);
236  H5Dclose(dataset_id2);
237  H5Dclose(dataset_id3);
238  H5Dclose(dataset_id4);
239 
240  H5Gclose(group_id);
241 
242  H5Fclose(file_id);
243 
244  if (CheMPS2::HAMILTONIAN_debugPrint) debugcheck();
245 
246 }
247 
248 void CheMPS2::Hamiltonian::CreateAndFillFromH5(const string file_parent, const string file_tmat, const string file_vmat){
249 
250  //The hdf5 file
251  hid_t file_id = H5Fopen(file_parent.c_str(), H5F_ACC_RDONLY, H5P_DEFAULT);
252 
253  //The data
254  hid_t group_id = H5Gopen(file_id, "/Data",H5P_DEFAULT);
255 
256  //The number of orbitals: Set L directly!
257  hid_t dataset_id1 = H5Dopen(group_id, "L", H5P_DEFAULT);
258  H5Dread(dataset_id1, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &L);
259 
260  //The group number: SymmInfo.setGroup()
261  hid_t dataset_id2 = H5Dopen(group_id, "nGroup", H5P_DEFAULT);
262  int nGroup_LOADH5;
263  H5Dread(dataset_id2, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &nGroup_LOADH5);
264  SymmInfo.setGroup( nGroup_LOADH5 );
265 
266  //The irrep of each orbital: Create and fill orb2irrep directly!
267  hid_t dataset_id3 = H5Dopen(group_id, "orb2irrep", H5P_DEFAULT);
268  orb2irrep = new int[ L ];
269  H5Dread(dataset_id3, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, orb2irrep);
270 
271  H5Dclose(dataset_id1);
272  H5Dclose(dataset_id2);
273  H5Dclose(dataset_id3);
274 
275  H5Gclose(group_id);
276 
277  H5Fclose(file_id);
278 
279  orb2indexSy = new int[ L ];
280  int nIrreps = SymmInfo.getNumberOfIrreps();
281  irrep2num_orb = new int[ nIrreps ];
282  for (int irrep = 0; irrep < nIrreps; irrep++){ irrep2num_orb[ irrep ] = 0; }
283  for (int orb = 0; orb < L; orb++){
284  orb2indexSy[ orb ] = irrep2num_orb[ orb2irrep[ orb ] ];
285  irrep2num_orb[ orb2irrep[ orb ] ]++;
286  }
287 
288  Tmat = new TwoIndex( SymmInfo.getGroupNumber(), irrep2num_orb );
289  Vmat = new FourIndex( SymmInfo.getGroupNumber(), irrep2num_orb );
290 
291  read(file_parent, file_tmat, file_vmat);
292 
293 }
294 
295 void CheMPS2::Hamiltonian::CreateAndFillFromFCIDUMP( const string fcidumpfile ){
296 
297  struct stat file_info;
298  const bool on_disk = (( fcidumpfile.length() > 0 ) && ( stat( fcidumpfile.c_str(), &file_info ) == 0 ));
299  if ( on_disk == false ){
300  cout << "CheMPS2::Hamiltonian : Unable to find FCIDUMP file " << fcidumpfile << "!" << endl;
301  }
302  assert( on_disk );
303 
304  cout << "CheMPS2::Hamiltonian : Reading FCIDUMP file " << fcidumpfile << endl;
305 
306  const int nIrreps = SymmInfo.getNumberOfIrreps();
307  int * psi2molpro = new int[ nIrreps ];
308  SymmInfo.symm_psi2molpro( psi2molpro );
309 
310  ifstream thefcidump( fcidumpfile.c_str() );
311  string line, part;
312  int pos, pos2;
313 
314  // Get the number of orbitals
315  getline( thefcidump, line ); // &FCI NORB= X,NELEC= Y,MS2= Z,
316  pos = line.find( "NORB" );
317  pos = line.find( "=", pos ); //1
318  pos2 = line.find( ",", pos ); //4
319  part = line.substr( pos+1, pos2-pos-1 );
320  L = atoi( part.c_str() );
321  if ( CheMPS2::HAMILTONIAN_debugPrint ){ cout << "The number of orbitals <<" << part << ">> or " << L << "." << endl; }
322 
323  // Get the orbital irreps in psi4 convention (XOR, see Irreps.h).
324  orb2irrep = new int[ L ];
325  getline( thefcidump, line ); // ORBSYM=A,B,C,D,
326  getline( thefcidump, part );
327  while ( part.find( "ISYM" ) == string::npos ){
328  pos = line.find("\n");
329  if ( pos != string::npos ){ line.erase( pos ); }
330  pos = part.find(" ");
331  if ( pos != string::npos ){ part.erase( pos, 1 ); }
332  line.append( part );
333  getline( thefcidump, part );
334  }
335  pos = line.find( "ORBSYM" );
336  pos = line.find( "=", pos ); //1
337  for ( int orb = 0; orb < L; orb++ ){
338  pos2 = line.find( ",", pos+1 ); //3
339  part = line.substr( pos+1, pos2-pos-1 );
340  const int molproirrep = atoi( part.c_str() );
341  if ( CheMPS2::HAMILTONIAN_debugPrint ){ cout << "This molpro irrep <<" << part << ">> or " << molproirrep << "." << endl; }
342  orb2irrep[ orb ] = -1;
343  for ( int irrep = 0; irrep < nIrreps; irrep++ ){
344  if ( molproirrep == psi2molpro[ irrep ] ){
345  orb2irrep[ orb ] = irrep;
346  }
347  }
348  assert( orb2irrep[ orb ] != -1 );
349  pos = pos2;
350  }
351 
352  getline( thefcidump, line ); // /
353  assert( line.size() < 16 );
354 
355  orb2indexSy = new int[ L ];
356  irrep2num_orb = new int[ nIrreps ];
357  for ( int cnt = 0; cnt < nIrreps; cnt++){ irrep2num_orb[cnt] = 0; }
358  for ( int cnt = 0; cnt < L; cnt++){
359  orb2indexSy[cnt] = irrep2num_orb[orb2irrep[cnt]];
360  irrep2num_orb[orb2irrep[cnt]]++;
361  }
362  Tmat = new TwoIndex( SymmInfo.getGroupNumber(), irrep2num_orb ); // Constructor ends with Clear(); call
363  Vmat = new FourIndex( SymmInfo.getGroupNumber(), irrep2num_orb ); // Constructor ends with Clear(); call
364 
365  // Read the Hamiltonian in
366  bool stop = false;
367  while ( stop == false ){
368 
369  getline( thefcidump, line ); // value i1 i2 i3 i4
370  pos = line.find( " " );
371  pos2 = line.find( "." );
372  pos2 = line.find( " ", pos2 );
373  part = line.substr( pos, pos2-pos );
374  if ( CheMPS2::HAMILTONIAN_debugPrint ){ cout << "Next line: <<" << part << ">> <<"; }
375  const double value = atof( part.c_str() );
376  pos = pos2;
377  while ( line.substr( pos, 1 ).compare(" ")==0 ){ pos++; }
378  pos2 = line.find( " ", pos );
379  part = line.substr( pos, pos2-pos );
380  if ( CheMPS2::HAMILTONIAN_debugPrint ){ cout << part << ">> <<"; }
381  const int index1 = atoi( part.c_str() );
382  pos = pos2;
383  while ( line.substr( pos, 1 ).compare(" ")==0 ){ pos++; }
384  pos2 = line.find( " ", pos );
385  part = line.substr( pos, pos2-pos );
386  if ( CheMPS2::HAMILTONIAN_debugPrint ){ cout << part << ">> <<"; }
387  const int index2 = atoi( part.c_str() );
388  pos = pos2;
389  while ( line.substr( pos, 1 ).compare(" ")==0 ){ pos++; }
390  pos2 = line.find( " ", pos );
391  part = line.substr( pos, pos2-pos );
392  if ( CheMPS2::HAMILTONIAN_debugPrint ){ cout << part << ">> <<"; }
393  const int index3 = atoi( part.c_str() );
394  pos = pos2;
395  while ( line.substr( pos, 1 ).compare(" ")==0 ){ pos++; }
396  part = line.substr( pos, line.size()-pos );
397  if ( CheMPS2::HAMILTONIAN_debugPrint ){ cout << part << ">>." << endl; }
398  const int index4 = atoi( part.c_str() );
399  if ( CheMPS2::HAMILTONIAN_debugPrint ){
400  cout << "Same line: " << value << " " << index1 << " " << index2 << " " << index3 << " " << index4 << endl;
401  }
402 
403  if ( index4 != 0 ){
404  setVmat( index1-1, index3-1, index2-1, index4-1, value ); // From chemists to physicist notation!
405  } else {
406  if ( index2 != 0 ){ setTmat( index1-1, index2-1, value ); }
407  else {
408  Econst = value;
409  stop = true;
410  }
411  }
412 
413  }
414 
415  if ( CheMPS2::HAMILTONIAN_debugPrint ){ debugcheck(); }
416 
417  delete [] psi2molpro;
418  thefcidump.close();
419 
420  cout << "CheMPS2::Hamiltonian : Finished reading FCIDUMP file " << fcidumpfile << endl;
421 
422 }
423 
424 void CheMPS2::Hamiltonian::readfock( const string fockfile, double * fockmx, const bool printinfo ) const{
425 
426 /************************
427  * FOCK file format *
428 ************************************
429  &FOCK NACT= X,
430  ORBSYM=X,X,X,X,
431  /
432  1.1234567890123456E-03 I J
433  1.1234567890123456E-03 I J
434 ************************************
435  Remark: ORBSYM are the MOLPRO irreps !!!
436 */
437 
438  struct stat file_info;
439  const bool on_disk = (( fockfile.length() > 0 ) && ( stat( fockfile.c_str(), &file_info ) == 0 ));
440  if ( on_disk == false ){
441  cout << "Unable to find the FOCK file " << fockfile << "!" << endl;
442  }
443  assert( on_disk );
444 
445  cout << "Reading FOCK file " << fockfile << endl;
446 
447  const int nIrreps = SymmInfo.getNumberOfIrreps();
448  int * psi2molpro = new int[ nIrreps ];
449  SymmInfo.symm_psi2molpro( psi2molpro );
450 
451  ifstream thedump( fockfile.c_str() );
452  string line, part;
453  int pos, pos2;
454 
455  // Double check the number of orbitals
456  getline( thedump, line ); // &FOCK NACT= X,
457  pos = line.find( "NACT" );
458  pos = line.find( "=", pos ); //1
459  pos2 = line.find( ",", pos ); //4
460  part = line.substr( pos+1, pos2-pos-1 );
461  const int LAS = atoi( part.c_str() );
462  if ( LAS != getL() ){
463  cout << "The number of orbitals in the FOCK file and Hamiltonian object does not match!" << endl;
464  assert( LAS == getL() );
465  }
466 
467  // Double check the orbital irreps
468  getline( thedump, line ); // ORBSYM=A,B,C,D,
469  getline( thedump, part );
470  while ( part.find( "/" ) == string::npos ){
471  pos = line.find( "\n" );
472  if ( pos != string::npos ){ line.erase( pos ); }
473  pos = part.find( " " );
474  if ( pos != string::npos ){ part.erase( pos, 1 ); }
475  line.append( part );
476  getline( thedump, part );
477  }
478  pos = line.find( "ORBSYM" );
479  pos = line.find( "=", pos ); //1
480  for ( int orb = 0; orb < LAS; orb++ ){
481  pos2 = line.find( ",", pos+1 ); //3
482  part = line.substr( pos+1, pos2-pos-1 );
483  const int molproirrep = atoi( part.c_str() );
484  if ( molproirrep != psi2molpro[ getOrbitalIrrep( orb ) ] ){
485  cout << "The irrep of orbital " << orb << " in the FOCK file and Hamiltonian object does not match!" << endl;
486  assert( molproirrep == psi2molpro[ getOrbitalIrrep( orb ) ] );
487  }
488  pos = pos2;
489  }
490 
491  // Read the FOCK matrix in
492  for ( int cnt = 0; cnt < LAS * LAS; cnt++ ){ fockmx[ cnt ] = 0.0; }
493  while ( getline( thedump, line ) ){ // value i1 i2
494 
495  if ( line.length() > 2 ){
496 
497  pos = line.find( " " );
498  pos2 = line.find( "." );
499  pos2 = line.find( " ", pos2 );
500  part = line.substr( pos, pos2-pos );
501  const double value = atof( part.c_str() );
502  pos = pos2;
503  while ( line.substr( pos, 1 ).compare(" ") == 0 ){ pos++; }
504  pos2 = line.find( " ", pos );
505  part = line.substr( pos, pos2-pos );
506  const int index1 = atoi( part.c_str() );
507  pos = pos2;
508  while ( line.substr( pos, 1 ).compare(" ") == 0 ){ pos++; }
509  pos2 = line.find( " ", pos );
510  part = line.substr( pos, pos2-pos );
511  const int index2 = atoi( part.c_str() );
512 
513  if ( printinfo ){
514  cout << " Processed FOCK( " << index1 << ", " << index2 << " ) = " << value << endl;
515  }
516 
517  const int orb1 = index1 - 1;
518  const int orb2 = index2 - 1;
519 
520  if ( getOrbitalIrrep( orb1 ) != getOrbitalIrrep( orb2 ) ){
521  cout << "In the FOCK file a specific value is given for orbitals " << index1 << " and " << index2 << " which have different irreps. This is not allowed!" << endl;
522  assert( getOrbitalIrrep( orb1 ) == getOrbitalIrrep( orb2 ) );
523  }
524 
525  fockmx[ orb1 + LAS * orb2 ] = value;
526  fockmx[ orb2 + LAS * orb1 ] = value;
527 
528  }
529  }
530 
531  delete [] psi2molpro;
532  thedump.close();
533 
534  cout << "Finished reading FOCK file " << fockfile << endl;
535 
536 }
537 
538 void CheMPS2::Hamiltonian::writeFCIDUMP( const string fcidumpfile, const int Nelec, const int TwoS, const int TargetIrrep ) const{
539 
540  int * psi2molpro = new int[ SymmInfo.getNumberOfIrreps() ];
541  SymmInfo.symm_psi2molpro( psi2molpro );
542 
543  FILE * capturing;
544  capturing = fopen( fcidumpfile.c_str(), "w" ); // "w" with fopen means truncate file
545  fprintf( capturing, " &FCI NORB= %d,NELEC= %d,MS2= %d,\n", getL(), Nelec, TwoS );
546  fprintf( capturing, " ORBSYM=" );
547  for (int orb=0; orb<getL(); orb++){
548  fprintf( capturing, "%d,", psi2molpro[getOrbitalIrrep(orb)] );
549  }
550  fprintf( capturing, "\n ISYM=%d,\n /\n", psi2molpro[TargetIrrep] );
551  delete [] psi2molpro;
552 
553  for (int p=0; p<getL(); p++){
554  for (int q=0; q<=p; q++){ // p>=q
555  const int irrep_pq = Irreps::directProd( getOrbitalIrrep(p), getOrbitalIrrep(q) );
556  for (int r=0; r<=p; r++){ // p>=r
557  for (int s=0; s<=r; s++){ // r>=s
558  const int irrep_rs = Irreps::directProd( getOrbitalIrrep(r), getOrbitalIrrep(s) );
559  if ( irrep_pq == irrep_rs ){
560  if ( ( p > r ) || ( ( p == r ) && ( q >= s ) ) ){
561  fprintf( capturing, " % 23.16E %3d %3d %3d %3d\n", getVmat(p,r,q,s), p+1, q+1, r+1, s+1 );
562  }
563  }
564  }
565  }
566  }
567  }
568 
569  for (int p=0; p<getL(); p++){
570  for (int q=0; q<=p; q++){ // p>=q
571  if ( getOrbitalIrrep(p) == getOrbitalIrrep(q) ){
572  fprintf( capturing, " % 23.16E %3d %3d %3d %3d\n", getTmat(p,q), p+1, q+1, 0, 0 );
573  }
574  }
575  }
576 
577  fprintf( capturing, " % 23.16E %3d %3d %3d %3d", getEconst(), 0, 0, 0, 0 );
578  fclose( capturing );
579  cout << "Created the FCIDUMP file " << fcidumpfile << "." << endl;
580 
581 }
582 
584 
585  cout << "Econst = " << Econst << endl;
586 
587  double test = 0.0;
588  double test2 = 0.0;
589  double test3 = 0.0;
590  for (int i=0; i<L; i++){
591  test3 += getTmat(i,i);
592  for (int j=0; j<L; j++){
593  test += getTmat(i,j);
594  if (i<=j) test2 += getTmat(i,j);
595  }
596  }
597  cout << "1-electron integrals: Trace : " << test3 << endl;
598  cout << "1-electron integrals: Sum over all elements : " << test << endl;
599  cout << "1-electron integrals: Sum over Tij with i<=j : " << test2 << endl;
600 
601  test = 0.0;
602  test2 = 0.0;
603  test3 = 0.0;
604  for (int i=0; i<L; i++){
605  test3 += getVmat(i,i,i,i);
606  for (int j=0; j<L; j++){
607  for (int k=0; k<L; k++){
608  for (int l=0; l<L; l++){
609  test += getVmat(i,j,k,l);
610  if ((i<=j) && (j<=k) && (k<=l)) test2 += getVmat(i,j,k,l);
611  }
612  }
613  }
614  }
615  cout << "2-electron integrals: Trace : " << test3 << endl;
616  cout << "2-electron integrals: Sum over all elements : " << test << endl;
617  cout << "2-electron integrals: Sum over Vijkl with i<=j<=k<=l : " << test2 << endl;
618 
619 }
620 
621 
void save(const std::string name) const
Save the FourIndex object.
Definition: FourIndex.cpp:314
bool setGroup(const int nGroup)
Set the group.
Definition: Irreps.cpp:50
void readfock(const string fockfile, double *fockmx, const bool printinfo) const
Read in a FOCK file.
int getNumberOfIrreps() const
Get the number of irreps for the currently activated group.
Definition: Irreps.cpp:94
void setTmat(const int index1, const int index2, const double val)
Set a Tmat element.
double get(const int irrep_i, const int irrep_j, const int irrep_k, const int irrep_l, const int i, const int j, const int k, const int l) const
Get an element.
Definition: FourIndex.cpp:175
int getOrbitalIrrep(const int nOrb) const
Get an orbital irrep number.
Definition: Hamiltonian.cpp:97
void debugcheck() const
Debug check certain elements and sums.
void setVmat(const int index1, const int index2, const int index3, const int index4, const double val)
Set a Vmat element.
void save(const string file_parent=HAMILTONIAN_ParentStorageName, const string file_tmat=HAMILTONIAN_TmatStorageName, const string file_vmat=HAMILTONIAN_VmatStorageName) const
Save the Hamiltonian.
void symm_psi2molpro(int *psi2molpro) const
Fill the array psi2molpro with the irrep conventions of molpro for the currently activated group...
Definition: Irreps.cpp:177
void setEconst(const double val)
Set the constant energy.
Definition: Hamiltonian.cpp:99
int getGroupNumber() const
Get the group number.
Definition: Irreps.cpp:66
virtual ~Hamiltonian()
Destructor.
Definition: Hamiltonian.cpp:83
void writeFCIDUMP(const string fcidumpfile, const int Nelec, const int TwoS, const int TargetIrrep) const
Write the Hamiltonian to a FCIDUMP file.
void read(const std::string name)
Load the TwoIndex object.
Definition: TwoIndex.cpp:140
void read(const string file_parent=HAMILTONIAN_ParentStorageName, const string file_tmat=HAMILTONIAN_TmatStorageName, const string file_vmat=HAMILTONIAN_VmatStorageName)
Load the Hamiltonian.
void set(const int irrep_i, const int irrep_j, const int irrep_k, const int irrep_l, const int i, const int j, const int k, const int l, const double val)
Set an element.
Definition: FourIndex.cpp:163
void addToVmat(const int index1, const int index2, const int index3, const int index4, const double val)
Add to Vmat element.
void set(const int irrep, const int i, const int j, const double val)
Set an element.
Definition: TwoIndex.cpp:64
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 get(const int irrep, const int i, const int j) const
Get an element.
Definition: TwoIndex.cpp:71
FourIndex * getVmat()
Get the pointer to the two-electron integrals.
int getNGroup() const
Get the group number.
Definition: Hamiltonian.cpp:95
double getEconst() const
Get the constant energy.
void save(const std::string name) const
Save the TwoIndex object.
Definition: TwoIndex.cpp:78
Hamiltonian(const int Norbitals, const int nGroup, const int *OrbIrreps)
Constructor.
Definition: Hamiltonian.cpp:39
const TwoIndex * getTmat()
Get the pointer to the one-electron integrals.
int getL() const
Get the number of orbitals.
Definition: Hamiltonian.cpp:93
void add(const int irrep_i, const int irrep_j, const int irrep_k, const int irrep_l, const int i, const int j, const int k, const int l, const double val)
Add a double to an element.
Definition: FourIndex.cpp:169
void read(const std::string name)
Load the FourIndex object.
Definition: FourIndex.cpp:372