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 |
None |
Optional dict with the damping parameters |
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 |
16.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!
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.0114416338147162
>>> atoms.calc.set(method="PBE")
{'method': 'PBE'}
>>> atoms.get_potential_energy()
-0.005358475432239303
>>> atoms.get_forces()
array([[-0. , -0. , 0.00296845],
[-0. , 0.00119152, -0.00148423],
[-0. , -0.00119152, -0.00148423]])