30 #include "FourIndex.h" 31 #include "Hamiltonian.h" 46 orb2irrep =
new int[L];
47 orb2indexSy =
new int[L];
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]]++;
67 SymmInfo.
setGroup( psi4groupnumber );
68 CreateAndFillFromFCIDUMP( filename );
75 CreateAndFillFromH5( main_file, file_tmat, file_vmat );
77 cout <<
"CheMPS2::Hamiltonian::Hamiltonian( false, const string , const string , const string ) was deprecated." << endl;
78 assert( fileh5 ==
true );
86 delete [] orb2indexSy;
87 delete [] irrep2num_orb;
105 assert( orb2irrep[index1]==orb2irrep[index2] );
106 Tmat->
set(orb2irrep[index1], orb2indexSy[index1], orb2indexSy[index2], val);
112 if (orb2irrep[index1]==orb2irrep[index2]){
113 return Tmat->
get(orb2irrep[index1], orb2indexSy[index1], orb2indexSy[index2]);
125 Vmat->
set(orb2irrep[index1], orb2irrep[index2], orb2irrep[index3], orb2irrep[index4], orb2indexSy[index1], orb2indexSy[index2], orb2indexSy[index3], orb2indexSy[index4], val);
132 Vmat->
add(orb2irrep[index1], orb2irrep[index2], orb2irrep[index3], orb2irrep[index4], orb2indexSy[index1], orb2indexSy[index2], orb2indexSy[index3], orb2indexSy[index4], val);
139 return Vmat->
get(orb2irrep[index1], orb2irrep[index2], orb2irrep[index3], orb2irrep[index4], orb2indexSy[index1], orb2indexSy[index2], orb2indexSy[index3], orb2indexSy[index4]);
150 Tmat->
save(file_tmat);
151 Vmat->
save(file_vmat);
154 hid_t file_id = H5Fcreate(file_parent.c_str(), H5F_ACC_TRUNC, H5P_DEFAULT, H5P_DEFAULT);
157 hid_t group_id = H5Gcreate(file_id,
"/Data", H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
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);
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);
170 H5Dwrite(dataset_id2, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &nGroup);
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);
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);
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);
201 Tmat->
read(file_tmat);
202 Vmat->
read(file_vmat);
205 hid_t file_id = H5Fopen(file_parent.c_str(), H5F_ACC_RDONLY, H5P_DEFAULT);
208 hid_t group_id = H5Gopen(file_id,
"/Data",H5P_DEFAULT);
211 hid_t dataset_id1 = H5Dopen(group_id,
"L", H5P_DEFAULT);
213 H5Dread(dataset_id1, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &Lagain);
217 hid_t dataset_id2 = H5Dopen(group_id,
"nGroup", H5P_DEFAULT);
219 H5Dread(dataset_id2, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &nGroup);
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] );
229 delete [] orb2irrepAgain;
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);
235 H5Dclose(dataset_id1);
236 H5Dclose(dataset_id2);
237 H5Dclose(dataset_id3);
238 H5Dclose(dataset_id4);
244 if (CheMPS2::HAMILTONIAN_debugPrint)
debugcheck();
248 void CheMPS2::Hamiltonian::CreateAndFillFromH5(
const string file_parent,
const string file_tmat,
const string file_vmat){
251 hid_t file_id = H5Fopen(file_parent.c_str(), H5F_ACC_RDONLY, H5P_DEFAULT);
254 hid_t group_id = H5Gopen(file_id,
"/Data",H5P_DEFAULT);
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);
261 hid_t dataset_id2 = H5Dopen(group_id,
"nGroup", H5P_DEFAULT);
263 H5Dread(dataset_id2, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &nGroup_LOADH5);
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);
271 H5Dclose(dataset_id1);
272 H5Dclose(dataset_id2);
273 H5Dclose(dataset_id3);
279 orb2indexSy =
new int[ L ];
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 ] ]++;
291 read(file_parent, file_tmat, file_vmat);
295 void CheMPS2::Hamiltonian::CreateAndFillFromFCIDUMP(
const string fcidumpfile ){
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;
304 cout <<
"CheMPS2::Hamiltonian : Reading FCIDUMP file " << fcidumpfile << endl;
307 int * psi2molpro =
new int[ nIrreps ];
310 ifstream thefcidump( fcidumpfile.c_str() );
315 getline( thefcidump, line );
316 pos = line.find(
"NORB" );
317 pos = line.find(
"=", pos );
318 pos2 = line.find(
",", pos );
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; }
324 orb2irrep =
new int[ L ];
325 getline( thefcidump, line );
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 ); }
333 getline( thefcidump, part );
335 pos = line.find(
"ORBSYM" );
336 pos = line.find(
"=", pos );
337 for (
int orb = 0; orb < L; orb++ ){
338 pos2 = line.find(
",", pos+1 );
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;
348 assert( orb2irrep[ orb ] != -1 );
352 getline( thefcidump, line );
353 assert( line.size() < 16 );
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]]++;
367 while ( stop ==
false ){
369 getline( thefcidump, line );
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() );
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() );
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() );
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() );
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;
404 setVmat( index1-1, index3-1, index2-1, index4-1, value );
406 if ( index2 != 0 ){
setTmat( index1-1, index2-1, value ); }
415 if ( CheMPS2::HAMILTONIAN_debugPrint ){
debugcheck(); }
417 delete [] psi2molpro;
420 cout <<
"CheMPS2::Hamiltonian : Finished reading FCIDUMP file " << fcidumpfile << endl;
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;
445 cout <<
"Reading FOCK file " << fockfile << endl;
448 int * psi2molpro =
new int[ nIrreps ];
451 ifstream thedump( fockfile.c_str() );
456 getline( thedump, line );
457 pos = line.find(
"NACT" );
458 pos = line.find(
"=", pos );
459 pos2 = line.find(
",", pos );
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() );
468 getline( thedump, line );
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 ); }
476 getline( thedump, part );
478 pos = line.find(
"ORBSYM" );
479 pos = line.find(
"=", pos );
480 for (
int orb = 0; orb < LAS; orb++ ){
481 pos2 = line.find(
",", pos+1 );
482 part = line.substr( pos+1, pos2-pos-1 );
483 const int molproirrep = atoi( part.c_str() );
485 cout <<
"The irrep of orbital " << orb <<
" in the FOCK file and Hamiltonian object does not match!" << endl;
492 for (
int cnt = 0; cnt < LAS * LAS; cnt++ ){ fockmx[ cnt ] = 0.0; }
493 while ( getline( thedump, line ) ){
495 if ( line.length() > 2 ){
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() );
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() );
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() );
514 cout <<
" Processed FOCK( " << index1 <<
", " << index2 <<
" ) = " << value << endl;
517 const int orb1 = index1 - 1;
518 const int orb2 = index2 - 1;
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;
525 fockmx[ orb1 + LAS * orb2 ] = value;
526 fockmx[ orb2 + LAS * orb1 ] = value;
531 delete [] psi2molpro;
534 cout <<
"Finished reading FOCK file " << fockfile << endl;
544 capturing = fopen( fcidumpfile.c_str(),
"w" );
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++){
550 fprintf( capturing,
"\n ISYM=%d,\n /\n", psi2molpro[TargetIrrep] );
551 delete [] psi2molpro;
553 for (
int p=0; p<
getL(); p++){
554 for (
int q=0; q<=p; q++){
556 for (
int r=0; r<=p; r++){
557 for (
int s=0; s<=r; 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 );
569 for (
int p=0; p<
getL(); p++){
570 for (
int q=0; q<=p; q++){
572 fprintf( capturing,
" % 23.16E %3d %3d %3d %3d\n",
getTmat(p,q), p+1, q+1, 0, 0 );
577 fprintf( capturing,
" % 23.16E %3d %3d %3d %3d",
getEconst(), 0, 0, 0, 0 );
579 cout <<
"Created the FCIDUMP file " << fcidumpfile <<
"." << endl;
585 cout <<
"Econst = " << Econst << endl;
590 for (
int i=0; i<L; i++){
592 for (
int j=0; j<L; j++){
594 if (i<=j) test2 +=
getTmat(i,j);
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;
604 for (
int i=0; i<L; 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++){
610 if ((i<=j) && (j<=k) && (k<=l)) test2 +=
getVmat(i,j,k,l);
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;
void save(const std::string name) const
Save the FourIndex object.
bool setGroup(const int nGroup)
Set the group.
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.
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.
int getOrbitalIrrep(const int nOrb) const
Get an orbital irrep number.
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...
void setEconst(const double val)
Set the constant energy.
int getGroupNumber() const
Get the group number.
virtual ~Hamiltonian()
Destructor.
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.
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.
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.
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's convent...
double get(const int irrep, const int i, const int j) const
Get an element.
FourIndex * getVmat()
Get the pointer to the two-electron integrals.
int getNGroup() const
Get the group number.
double getEconst() const
Get the constant energy.
void save(const std::string name) const
Save the TwoIndex object.
Hamiltonian(const int Norbitals, const int nGroup, const int *OrbIrreps)
Constructor.
const TwoIndex * getTmat()
Get the pointer to the one-electron integrals.
int getL() const
Get the number of orbitals.
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.
void read(const std::string name)
Load the FourIndex object.