vdW-DF and BEEF-vdW

GPAW supports vdW-DF functionals through a built-in interface as well as through the external libvdwxc library. Note that these use different kernels and hence will yield slightly different results.

Several vdW-DF [1] type XC functionals are implemented self-consistently in GPAW, and also the BEEF-vdW [2] density functional. The vdW-DF variants include vdW-DF [1], [3], vdW-DF2 [4], vdW-DF-cx, [5], optPBE-vdW [6], optB88-vdW [6], and C09-vdW [7].

Of these, vdW-DF-cx is available only through libvdwxc. The spin-polarized generalization of the vdW-DF functionals, [9], is also only available with libvdwxc.

The self-consistent implementation uses the Perez-Soler [8] FFT algorithm to evaluate the total energy and potential of the Rutgers-Chalmers nonlocal correlation, which is originally a six dimensional integral in real space. However, a non self-consistent method which directly sums up the real-space integral is also available.

Doing a vdW-DF calculation

The self-consistent FFT method is highly recommended over the real-space method. Often, the vdW-DF electron density will be very similar to an ordinary GGA density, so non self-consistent evaluations of a vdW-DF type total energy using the FFT method is often ok. However, vdW-DF forces obviously require a self-consistent potential.

As the examples below illustrate, FFT-based vdW-DF calculations are most easily done by setting e.g. “xc=’vdW-DF’” in the GPAW calculator object. However, parameters of the FFT algorithm can be assigned non-default values by importing the vdW-DF base class.

For larger systems, the van der Waals functionals may be computationally expensive. Consider using libvdwxc which typically increases the efficiency of the van der Waals evaluation by an order of magnitude, and parallelizes to any desired system size.

Selfconsistent vdW-DF calculations

>>> from ase import *
>>> from gpaw import GPAW
>>> vdw = 'vdW-DF'
>>> atoms = ...
>>> calc = GPAW(xc=vdw, ...)
>>> atoms.calc = calc
>>> e = atoms.get_potential_energy()

Perturbative vdW-DF calculations (non self-consistent)

>>> from gpaw import GPAW
>>> xc = 'vdW-DF'
>>> calc = GPAW('input.gpw')
>>> GGA_energy = calc.get_potential_energy()
>>> vdWDF_diff = calc.get_xc_difference(xc)
>>> vdWDF_energy = GGA_energy + vdWDF_diff

In the above examples, other vdW-DF type functionals can be used by substituting ‘vdW-DF2’, ‘vdW-DF-cx’ (if GPAW is compiled with libvdwxc), ‘optPBE-vdW’, ‘optB88-vdW’, or ‘C09-vdW’ for ‘vdW-DF’. To explicitly use the faster libvdwxc backend, use e.g. xc={'name': 'vdW-DF', 'backend': 'libvdwxc'}. t libvdwxc uses a different kernel parametrization, which will slightly affect calculated values.

Non-default FFT parameters for vdW-DF calculations

A number of parameters determine the spline interpolation of the vdW-DF nonlocal kernel. These may be assigned non-default values if the vdW-DF base class is explicitly initialized with new settings. The example below redefines the number of interpolating cubic splines (Nalpha) used in a vdW-DF2 calculation.

>>> from ase import *
>>> from gpaw import GPAW
>>> from gpaw.xc.vdw import VDWFunctional
>>> vdw = VDWFunctional('vdW-DF2', Nalpha=24)
>>> atoms = ...
>>> calc = GPAW(xc=vdw, ...)
>>> atoms.calc = calc
>>> e = atoms.get_potential_energy()

Real-space method vdW-DF

It is also possible to use the much slower real-space method for non self-consistent evaluations of the nonlocal correlation energy, which might make sense for (very) small systems. To use the real-space method one must import a class and set a few parameters:

>>> from gpaw.xc.vdw import VDWFunctional
>>> vdw = VDWFunctional('vdW-DF', fft=False, nspins=1, ncut=0.0005)

where nspins=1 is for spin-paired systems and nspins=2 is used for spin-polarized calculations. A cutoff, ncut, defines how small a density must be in order not to be included in the 6D integral.

BEEF-vdW functional

The BEEF-vdW density functional uses the vdW-DF2 nonlocal correlation energy and potential. It is implemented selfconistently in GPAW. Furthermore, the BEEF-vdW constructions allows the user to calculate an estimate of the error to be expected on the quantity calculated self-consistently with BEEF-vdW (i.e. an error estimate on relative energies, not on total energies). This estimate stems from non self-consistently applying an ensemble of XC functionals to BEEF-vdW electron densities. The ensemble error estimate is then computed from the variance of the ensemble predictions of the quantity of interest.

Below is an example which calculates the BEEF-vdW binding energy of molecular H2 (E_bind), as well as an ensemble estimate of the binding energy error (dE_bind)

>>> from ase import *
>>> from gpaw import GPAW
>>> from ase.dft.bee import BEEFEnsemble
>>> xc = 'BEEF-vdW'
>>> h2 = Atoms('H2',[[0.,0.,0.],[0.,0.,0.75]])
>>> h2.center(vacuum=3)
>>> cell = h2.get_cell()
>>> calc = GPAW(mode='fd', xc=xc)
>>> h2.calc = calc
>>> e_h2 = h2.get_potential_energy()
>>> ens = BEEFEnsemble(calc)
>>> de_h2 = ens.get_ensemble_energies()
>>> del h2, calc, ens
>>> h = Atoms('H')
>>> h.set_cell(cell)
>>> h.center()
>>> calc = GPAW(mode='fd', xc=xc)
>>> h.calc = calc
>>> e_h = h.get_potential_energy()
>>> ens = BEEFEnsemble(calc)
>>> de_h = ens.get_ensemble_energies()
>>> E_bind = 2*e_h - e_h2
>>> dE_bind = 2*de_h[:] - de_h2[:]
>>> dE_bind = dE_bind.std()

Note that the BEEFEnsemble module has recently been moved from GPAW to the ASE package. The default number of ensemble XC functionals is 2000, for which well-converged error estimates should be ensured. Therefore, “de_h2” and “de_h” in the example are both arrays of 2000 perturbations of a BEEF-vdW total energy. The syntax “ens.get_ensemble_energies(N)” changes this number to N. The calculator object input to the BEEFEnsemble class could of course stem from a restarted GPAW calculation.

It is very important to calculate the ensemble statistics correctly. Computing the standard deviation of each array of total energy perturbations makes little sense, only the standard deviation of the relative energy perturbations should be used for the BEEF-vdW ensemble error estimates on a quantity.