ASE Support#
ASE calculator implementation
for the s-dftd3 program.
This module provides a basic single point calculator implementations
to integrate the s-dftd3 API into existing ASE workflows.
To use DFTD3 as dispersion correction the ase.calculators.mixing
module can be used to combine DFTD3 with a DFT calculator using
the SumCalculator.
Supported properties by this calculator are:
energy (free_energy)
forces
stress
Supported keywords are
Keyword |
Default |
Description |
|---|---|---|
method |
None |
Method to calculate dispersion for |
damping |
None |
Damping function to use |
params_tweaks |
{} |
Optional dict with the damping parameters |
realspace_cutoff |
{} |
Optional dict to override cutoff values |
ghost_atoms |
None |
Disable dispersion contributions from these atoms |
cache_api |
True |
Reuse generate API objects (recommended) |
The params_tweaks dict contains the damping parameters, at least s8, a1 and a2 must be provided
Tweakable parameter |
Default |
Description |
|---|---|---|
s6 |
1.0 |
Scaling of the dipole-dipole dispersion |
s8 |
None |
Scaling of the dipole-quadrupole dispersion |
s9 |
1.0 |
Scaling of the three-body dispersion energy |
a1 |
None |
Scaling of the critical radii |
a2 |
None |
Offset of the critical radii |
alp |
14.0 |
Exponent of the zero damping (ATM only) |
Either method or s8, a1 and a2 must be provided, s9 can be used to overwrite the ATM scaling if the method is provided in the model. Disabling the three-body dispersion (s9=0.0) changes the internal selection rules for damping parameters of a given method and prefers special two-body only damping parameters if available!
The realspace cutoff parameters allow adjusting the distance values for which interactions are considered
Realspace cutoff |
Default |
Description |
|---|---|---|
disp2 |
60 * Bohr |
Pairwise dispersion interactions |
disp3 |
40 * Bohr |
Triple dispersion interactions |
cn |
40 * Bohr |
Coordination number counting |
width2 |
0 * Bohr |
Smooth cutoff width for pairwise dispersion |
width3 |
0 * Bohr |
Smooth cutoff width for triple dispersion |
Values provided in the dict are expected to be in Angstrom. When providing values in Bohr multiply the inputs by the ase.units.Bohr constant.
Example
>>> from ase.build import molecule
>>> from dftd3.ase import DFTD3
>>> atoms = molecule('H2O')
>>> atoms.calc = DFTD3(method="TPSS", damping="d3bj")
>>> atoms.get_potential_energy()
-0.011441640787279111
>>> atoms.calc.set(method="PBE")
{'method': 'PBE'}
>>> atoms.get_potential_energy()
-0.009781920198843976
>>> atoms.get_forces()
array([[-0. , -0. , 0.00009572],
[-0. , -0.00004060, -0.00004786],
[-0. , 0.00004060, -0.00004786]])