Skip to content

Commit

Permalink
use solvent model for unrestricted
Browse files Browse the repository at this point in the history
  • Loading branch information
wxj6000 committed Feb 6, 2024
1 parent eaf83a4 commit f64fcf3
Show file tree
Hide file tree
Showing 4 changed files with 21 additions and 3 deletions.
5 changes: 3 additions & 2 deletions gpu4pyscf/solvent/pcm.py
Original file line number Diff line number Diff line change
Expand Up @@ -42,7 +42,7 @@ def pcm_for_scf(mf, solvent_obj=None, dm=None):
# Inject PCM to SCF, TODO: add it to other methods later
from gpu4pyscf import scf
scf.hf.RHF.PCM = pcm_for_scf

scf.uhf.UHF.PCM = pcm_for_scf
# TABLE II, J. Chem. Phys. 122, 194110 (2005)
XI = {
6: 4.84566077868,
Expand Down Expand Up @@ -313,7 +313,8 @@ def _get_vind(self, dms):

nao = dms.shape[-1]
dms = dms.reshape(-1,nao,nao)

if dms.shape[0] == 2:
dms = (dms[0] + dms[1]).reshape(-1,nao,nao)
K = self._intermediates['K']
R = self._intermediates['R']
v_grids_e = self._get_v(dms)
Expand Down
2 changes: 1 addition & 1 deletion gpu4pyscf/solvent/smd.py
Original file line number Diff line number Diff line change
Expand Up @@ -36,7 +36,7 @@ def smd_for_scf(mf, solvent_obj=None, dm=None):
# Inject PCM to SCF, TODO: add it to other methods later
from gpu4pyscf import scf
scf.hf.RHF.SMD = smd_for_scf

scf.uhf.UHF.SMD = smd_for_scf
hartree2kcal = 627.5
# see https://pubs.acs.org/doi/epdf/10.1021/jp810292n

Expand Down
10 changes: 10 additions & 0 deletions gpu4pyscf/solvent/tests/test_pcm.py
Original file line number Diff line number Diff line change
Expand Up @@ -83,6 +83,16 @@ def test_SSVPE(self):
print(f"Energy error in SS(V)PE: {numpy.abs(e_tot - -74.9689577454)}")
assert numpy.abs(e_tot - -74.9689577454) < 1e-9

def test_uhf(self):
cm = pcm.PCM(mol)
cm.eps = epsilon
cm.verbose = 0
cm.lebedev_order = 29
cm.method = 'IEF-PCM'
mf = scf.UHF(mol).PCM(cm)
e_tot = mf.kernel()
assert numpy.abs(e_tot - -74.96901113434953) < 1e-9

if __name__ == "__main__":
print("Full Tests for PCMs")
unittest.main()
7 changes: 7 additions & 0 deletions gpu4pyscf/solvent/tests/test_smd.py
Original file line number Diff line number Diff line change
Expand Up @@ -68,6 +68,13 @@ def test_smd_water(self):
e_tot = mf.kernel()
assert numpy.abs(e_tot - -76.0756052903) < 2e-4

def test_uhf(self):
mf = scf.UHF(mol)
mf = mf.SMD()
mf.with_solvent.solvent = 'water'
mf.with_solvent.sasa_ng = 590
e_tot = mf.kernel()
assert numpy.abs(e_tot - -76.07550951172617) < 2e-4
# TODO: add more test for other molecules

if __name__ == "__main__":
Expand Down

0 comments on commit f64fcf3

Please sign in to comment.