Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

wB97m-D3BJ #95

Open
sef43 opened this issue Jan 31, 2024 · 6 comments
Open

wB97m-D3BJ #95

sef43 opened this issue Jan 31, 2024 · 6 comments

Comments

@sef43
Copy link

sef43 commented Jan 31, 2024

Hello,

is it possible to use the functional wB97m-D3BJ in gpu4pyscf?

It is wB97m-V with the VV10 replaced with D3(BJ). https://pubs.acs.org/doi/abs/10.1021/acs.jctc.8b00842

It is available in PSI4 but I would like to use PySCF for the GPU acceleration

@wxj6000
Copy link
Collaborator

wxj6000 commented Feb 1, 2024

Yes, it is possible. There will some script to work around the xc functional names. In libXC, there is only one functional for wb97m-v, which includes VV10. But the name wb97m-v is not available in dftd3. Here is the script to replace VV10 with D3(BJ).

import numpy as np
import cupy
import pyscf
from pyscf import lib
from gpu4pyscf.dft import rks, libxc
lib.num_threads(8)

atom = '''
O       0.0000000000    -0.0000000000     0.1174000000
H      -0.7570000000    -0.0000000000    -0.4696000000
H       0.7570000000     0.0000000000    -0.4696000000
'''

xc = 'wb97m-v'
bas = 'def2-tzvpp'
scf_tol = 1e-10
max_scf_cycles = 50
grids_level = 5

mol = pyscf.M(atom=atom, basis=bas)

mol.verbose = 3
mf = rks.RKS(mol, xc=xc).density_fit()

# turn off nonlocal correction
mf.nlc = None
mf._numint.libxc.is_nlc = lambda x: False
mf.grids.level = grids_level
mf.grids.atom_grid = (99,590)
mf.conv_tol = scf_tol
mf.max_cycle = max_scf_cycles
e_dft = mf.kernel()

# calculate dispersion correction
from gpu4pyscf.lib import dftd3
dftd3_model = dftd3.DFTD3Dispersion(mol, xc='wb97m', version='d3bj')
e_disp = res = dftd3_model.get_dispersion()['energy']
e_tot = e_dft + e_disp

print(f"total energy = {e_tot}")

@sef43
Copy link
Author

sef43 commented Feb 1, 2024

thank you very much for the clear example! I will try it out

@sef43
Copy link
Author

sef43 commented Jun 28, 2024

Hello. I only now got around to trying this. I have used the new syntax I found in the example: https://github.com/pyscf/gpu4pyscf/blob/master/examples/01-h2o_with_dispersion.py

I get an error when computing the gradient. How do I also compute the gradient with the D3bj?

here is the error:

Traceback (most recent call last):
  File "run_test.py", line 36, in <module>
    g_dft = g.kernel()
            ^^^^^^^^^^
  File "/scratch/users/sfarr/miniconda3/envs/gpu4pyscf_new/lib/python3.11/site-packages/pyscf/grad/rhf.py", line 415, in kernel
    de = self.grad_elec(mo_energy, mo_coeff, mo_occ, atmlst)
         ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
  File "/scratch/users/sfarr/miniconda3/envs/gpu4pyscf_new/lib/python3.11/site-packages/gpu4pyscf/grad/rhf.py", line 523, in grad_elec
    dvhf = mf_grad.get_veff(mol, dm0)
           ^^^^^^^^^^^^^^^^^^^^^^^^^^
  File "/scratch/users/sfarr/miniconda3/envs/gpu4pyscf_new/lib/python3.11/site-packages/gpu4pyscf/df/grad/rks.py", line 55, in get_veff
    raise NotImplementedError
NotImplementedError

and here is my script:

import numpy as np
import pyscf
from gpu4pyscf.dft import rks


atom ='''
O       0.0000000000    -0.0000000000     0.1174000000
H      -0.7570000000    -0.0000000000    -0.4696000000
H       0.7570000000     0.0000000000    -0.4696000000
'''

mol = pyscf.M(
    atom=atom,                         # water molecule
    basis='def2-tzvppd',                # basis set
    verbose=1                          # control the level of print info
    )

mf_GPU = rks.RKS(                      # restricted Kohn-Sham DFT
    mol,                               # pyscf.gto.object
    xc='wB97m-V'                         # xc funtionals, such as pbe0, wb97m-v, tpss,
    ).density_fit()                    # density fitting


# turn off nlc
mf_GPU.nlc = None
mf_GPU._numint.libxc.is_nlc = lambda x: False

mf_GPU.disp = 'd3bj'

# Compute Energy
e_dft = mf_GPU.kernel()
print(f"total energy = {e_dft}")

# Compute Gradient
g = mf_GPU.nuc_grad_method()
g_dft = g.kernel()
print(g_dft)

Thank you

@sef43
Copy link
Author

sef43 commented Jul 1, 2024

can now be done with:

mf_GPU = rks.RKS(                     
    mol,                               
    xc='wB97m-d3bj'                        
    ).density_fit() 

thank you @wxj6000

@sef43 sef43 closed this as completed Jul 1, 2024
@wxj6000 wxj6000 reopened this Jul 1, 2024
@wxj6000
Copy link
Collaborator

wxj6000 commented Jul 1, 2024

@sef43 Yes, PySCF v2.6 and GPU4PySCF v0.8 support RKS(mol, xc='wb97m-d3bj'). I deleted the comment before because I am planning to do more tests. I will add more unit tests for the new features.

@wxj6000
Copy link
Collaborator

wxj6000 commented Jul 3, 2024

#173
The PR mentioned above.

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
None yet
Projects
None yet
Development

No branches or pull requests

2 participants