TS09 and D4 corrections with ASE

TS09 and D4 are atomic-charge dependent dispersion corrections (see TS09 PRL paper and D4 homepage for the refs). The D4 code is available at github. According to GPAW documentation, TS09 and D4 show for the S26 test set smaller mean deviation than vdW-DF. Herewith, D4 correction does not depend on the actual calculation as it is added to the calculated energy.

Here is how D4 correction can be added with ASE (see Readme) after installing it (for example, as conda install -c conda-forge dftd4 dftd4-python):

from ase.build import molecule 
from ase.calculators.mixing import SumCalculator 
from ase.optimize import BFGS
from dftd4.ase import DFTD4 
from gpaw import GPAW 

atoms = molecule('H2O') 
atoms.center(vacuum=4)

gpaw = GPAW(txt='H2O_D4.txt',xc='PBE') 
atoms.calc = SumCalculator([DFTD4(method='PBE'), gpaw])

#atoms.get_potential_energy()
opt = BFGS(atoms,trajectory='H2O_D4.traj', logfile='H2O_D4.log')
opt.run(fmax=0.05)

Let me stress that before choosing TS09 or D4 one should consider all pro and contra. TS09 method used Hirshfeld charges while D4 uses the electronegativity equilibration method to obtain charges. The former naturally accounts for the interfacial charge transfer while the latter does not. The TS09 correction requires vdW radii and is implemented for a limited set on functionals (see ASE code), like PBE, RPBE, and BLYP. The D4 correction supports much more functionals (see parameters). Regarding the vdW radii values for TS09 bare in mind that there are four data sources – one in GPAW, two in ASE and one more in ASE.

Here is how TS09 correction can be added with ASE and GPAW:

from ase.build import molecule
from ase.calculators.vdwcorrection import vdWTkatchenko09prl
from ase.data.vdw_alvarez import vdw_radii
from ase.optimize import BFGS
from gpaw.analyse.hirshfeld import HirshfeldPartitioning
from gpaw.analyse.vdwradii import vdWradii
from gpaw import GPAW

atoms = molecule('H2O')
atoms.center(vacuum=4)

gpaw = GPAW(txt='H2O_TS.txt',xc='PBE')
atoms.calc = vdWTkatchenko09prl(HirshfeldPartitioning(gpaw), vdWradii(atoms.get_chemical_symbols(), 'PBE'))

#atoms.get_potential_energy()
opt = BFGS(atoms,trajectory='H2O_TS.traj', logfile='H2O_TS.log')
opt.run(fmax=0.0)

P.S. Note that the TS09 and D4 energies are no outputted to the H2O.txt. They are written to the log-file.

Bader and Hirshfeld charges with python

Bader analysis is a fast and simple way of getting atomic charges. It is especially useful for periodic calculations. The analysis can be done on the fly using pybader tool from pypi.org/project/pybader. I recommend using it within conda environments.

The installation is straightforward:

pip install pybader

The usage is less obvious. Here is a function for obtaining xmol-type xyz that can obtained with GPAW and visualized with ASE:


def xyzb(atoms, filename, nCPU):
  from pybader.io import gpaw
  from pybader.interface import Bader
  import os
  bader = Bader(*gpaw.read_obj(atoms.calc)) # read ASE object 'atoms'
  bader(threads=nCPU)                       # specify the number of CPUs
  f = open('{0}.xyz'.format(filename), 'w') # set xmol format
  b = bader.atoms_charge                    # get number of electrons per atom
  n = atoms.get_atomic_numbers()            # get atomic numbers
  a = atoms.get_chemical_symbols()          # get chemical symbols
  p = atoms.get_positions()                 # get positions of the atoms
  f.write('{0}\n'.format(len(a)))
  f.write('Properties=species:S:1:pos:R:3:charge:R:1\n') # ensure compatibility with ASE
  for i in range(len(a)):                   # print symbol, positions, and charge
    s = '{0}'.format(a[i])
    x = '{0:.6f}'.format(round(p[i][0],6))
    y = '{0:.6f}'.format(round(p[i][1],6))
    z = '{0:.6f}'.format(round(p[i][2],6))
    c = '{0:.3f}'.format(round(float(n[i]) - float(b[i]),3))
    f.write('{0:<4}{1:>16}{2:>16}{3:>16}{4:>10}\n'.format(s,x,y,z,c))
  f.close()
  del bader
  os.remove('bader.p')

Similarly one can obtain xmol-type xyz file with Hirshfeld charges:


def xyzh(atoms, filename):
  from gpaw.analyse.hirshfeld import HirshfeldPartitioning
  f = open('{0}.xyz'.format(filename), 'w') # set xmol format
  a = atoms.get_chemical_symbols()          # get chemical symbols
  p = atoms.get_positions()                 # get positions of the atoms
  h = HirshfeldPartitioning(atoms.calc).get_charges()
  f.write('{0}\n'.format(len(a)))
  f.write('Properties=species:S:1:pos:R:3:charge:R:1\n') # ensure compatibility with ASE
  for i in range(len(a)):                   # print symbol, positions, and charge
    s = '{0}'.format(a[i])
    x = '{0:.6f}'.format(round(p[i][0],6))
    y = '{0:.6f}'.format(round(p[i][1],6))
    z = '{0:.6f}'.format(round(p[i][2],6))
    c = '{0:.3f}'.format(round(h[i],3))
    f.write('{0:<4}{1:>16}{2:>16}{3:>16}{4:>10}\n'.format(s,x,y,z,c))
  f.close()

In the LCAO mode of GPAW one can also get the Mulliken charges. Test before using:


def xyzm(atoms, filename):
  from gpaw.lcao.tools import get_mulliken
  f = open('{0}.xyz'.format(filename), 'w') # set xmol format
  a = atoms.get_chemical_symbols()          # get chemical symbols
  p = atoms.get_positions()                 # get positions of the atoms
  m = get_mulliken(atoms.calc, range(len(a)))
  f.write('{0}\n'.format(len(a)))
  f.write('Properties=species:S:1:pos:R:3:charge:R:1\n') # ensure compatibility with ASE
  for i in range(len(a)):                   # print symbol, positions, and charge
    s = '{0}'.format(a[i])
    x = '{0:.6f}'.format(round(p[i][0],6))
    y = '{0:.6f}'.format(round(p[i][1],6))
    z = '{0:.6f}'.format(round(p[i][2],6))
    c = '{0:.3f}'.format(round(m[i],3))
    f.write('{0:<4}{1:>16}{2:>16}{3:>16}{4:>10}\n'.format(s,x,y,z,c))
  f.close()