Skip to content

Commit

Permalink
Fix chkfile initial guess problem for ROHF and ROKS (issue pyscf#1986)
Browse files Browse the repository at this point in the history
  • Loading branch information
sunqm committed Jan 1, 2024
1 parent 0bc590f commit 5aef8c3
Show file tree
Hide file tree
Showing 3 changed files with 58 additions and 5 deletions.
16 changes: 16 additions & 0 deletions pyscf/dft/test/test_he.py
Original file line number Diff line number Diff line change
Expand Up @@ -14,6 +14,7 @@
# limitations under the License.

import unittest
import tempfile
import numpy
from pyscf import gto
from pyscf import lib
Expand Down Expand Up @@ -191,6 +192,21 @@ def test_convert(self):
self.assertTrue(isinstance(udks.to_dhf(), scf.dhf.DHF))
self.assertTrue(isinstance(udks.to_dks('pbe'), dft.dks.DKS))

# issue 1986
def test_init_guess_chkfile(self):
with tempfile.NamedTemporaryFile() as tmpf:
mol = gto.M(atom='He 0 0 0', basis='631g', charge=1, spin=1)
mf = dft.RKS(mol)
mf.chkfile = tmpf.name
e1 = mf.kernel()
mf = dft.RKS(mol)
mf.init_guess = 'chkfile'
mf.chkfile = tmpf.name
mf.max_cycle = 1
e2 = mf.kernel()
self.assertAlmostEqual(e1, e2, 9)


if __name__ == "__main__":
print("Full Tests for He")
unittest.main()
24 changes: 20 additions & 4 deletions pyscf/scf/hf.py
Original file line number Diff line number Diff line change
Expand Up @@ -656,12 +656,28 @@ def init_guess_by_chkfile(mol, chkfile_name, project=None):
'''Read the HF results from checkpoint file, then project it to the
basis defined by ``mol``
Kwargs:
project : None or bool
Whether to project chkfile's orbitals to the new basis. Note when
the geometry of the chkfile and the given molecule are very
different, this projection can produce very poor initial guess.
In PES scanning, it is recommended to switch off project.
If project is set to None, the projection is only applied when the
basis sets of the chkfile's molecule are different to the basis
sets of the given molecule (regardless whether the geometry of
the two molecules are different). Note the basis sets are
considered to be different if the two molecules are derived from
the same molecule with different ordering of atoms.
Returns:
Density matrix, 2D ndarray
'''
from pyscf.scf import uhf
dm = uhf.init_guess_by_chkfile(mol, chkfile_name, project)
return dm[0] + dm[1]
mo_coeff = dm.mo_coeff[0]
mo_occ = dm.mo_occ[0] + dm.mo_occ[1]
return lib.tag_array(dm[0] + dm[1], mo_coeff=mo_coeff, mo_occ=mo_occ)


def get_init_guess(mol, key='minao'):
Expand Down Expand Up @@ -1358,11 +1374,9 @@ def __call__(self, mol_or_geom, **kwargs):
dm0 = kwargs.pop('dm0')
elif self.mo_coeff is None:
dm0 = None
elif self.chkfile and h5py.is_hdf5(self.chkfile):
dm0 = self.from_chk(self.chkfile)
else:
dm0 = None
# dm0 form last calculation cannot be used in the current
# dm0 form last calculation may not be used in the current
# calculation if a completely different system is given.
# Obviously, the systems are very different if the number of
# basis functions are different.
Expand All @@ -1371,6 +1385,8 @@ def __call__(self, mol_or_geom, **kwargs):
# last calculation.
if numpy.array_equal(self._last_mol_fp, mol.ao_loc):
dm0 = self.make_rdm1()
elif self.chkfile and h5py.is_hdf5(self.chkfile):
dm0 = self.from_chk(self.chkfile)
self.mo_coeff = None # To avoid last mo_coeff being used by SOSCF
e_tot = self.kernel(dm0=dm0, **kwargs)
self._last_mol_fp = mol.ao_loc
Expand Down
23 changes: 22 additions & 1 deletion pyscf/scf/rohf.py
Original file line number Diff line number Diff line change
Expand Up @@ -49,7 +49,28 @@ def init_guess_by_atom(mol):

init_guess_by_huckel = uhf.init_guess_by_huckel
init_guess_by_mod_huckel = uhf.init_guess_by_mod_huckel
init_guess_by_chkfile = uhf.init_guess_by_chkfile

def init_guess_by_chkfile(mol, chkfile_name, project=None):
'''Read SCF chkfile and make the density matrix for ROHF initial guess.
Kwargs:
project : None or bool
Whether to project chkfile's orbitals to the new basis. Note when
the geometry of the chkfile and the given molecule are very
different, this projection can produce very poor initial guess.
In PES scanning, it is recommended to switch off project.
If project is set to None, the projection is only applied when the
basis sets of the chkfile's molecule are different to the basis
sets of the given molecule (regardless whether the geometry of
the two molecules are different). Note the basis sets are
considered to be different if the two molecules are derived from
the same molecule with different ordering of atoms.
'''
dm = uhf.init_guess_by_chkfile(mol, chkfile_name, project)
mo_coeff = dm.mo_coeff[0]
mo_occ = dm.mo_occ[0] + dm.mo_occ[1]
return lib.tag_array(dm, mo_coeff=mo_coeff, mo_occ=mo_occ)

def get_fock(mf, h1e=None, s1e=None, vhf=None, dm=None, cycle=-1, diis=None,
diis_start_cycle=None, level_shift_factor=None, damp_factor=None):
Expand Down

0 comments on commit 5aef8c3

Please sign in to comment.