Usage in ASE#

The Python API of tblite natively support integration with the atomic simulation environment (ASE). By constructing a calculator most functionality of ASE is readily available. For details on building the Python API checkout the installation guide.

Example

The calculator can be used like other ASE calculators by attaching it to an atoms object.

>>> from ase.build import molecule
>>> from tblite.ase import TBLite
>>> atoms = molecule("H2O")
>>> atoms.calc = TBLite(method="GFN2-xTB")
>>> atoms.get_potential_energy()  # in eV
-137.96777594361677

The calculator will perform single point calculations as necessary, if properties for the same atoms are requested the properties will be returned without additional single point calculations.

>>> atoms.get_dipole_moment()  # in Debye
array([-7.83154860e-17,  1.10157136e-16, -4.76020974e-01])
>>> atoms.get_forces()  # in eV/A
array([[-3.31859412e-15,  7.83871951e-15, -7.49527168e-01],
       [ 3.08351626e-15, -1.53506496e-01,  3.74763584e-01],
       [ 2.35077857e-16,  1.53506496e-01,  3.74763584e-01]])

Changed properties of the atoms object will lead to an reevalution of the single point, for example updating the initial charges of the atoms to set a total charge.

>>> atoms.set_initial_charges([0, 0, -1])
>>> atoms.get_potential_energy()  # in eV
-131.50015843891816

Some properties, like the total charge, can also be set via the parameters of the calculator.

>>> atoms = molecule("H2O")
>>> atoms.calc = TBLite(method="GFN2-xTB", charge=-1)
>>> atoms.get_potential_energy()  # in eV
-131.50015843891822

Finally, the parameters of the calculator can be updated at any time using the set method, which resets the calculation as necessary.

>>> atoms.calc.set(method="GFN1-xTB", charge=0)
>>> atoms.get_potential_energy()  # in eV
-156.96750590276173
>>> atoms.get_dipole_moment()  # in Debye
array([ 3.71500789e-16,  3.08439980e-16, -5.97048616e-01])
>>> atoms.get_forces()  # in eV/A
array([[ 2.88997963e-15,  5.70620869e-15, -7.84252969e-01],
       [-4.90812472e-15, -2.28835628e-01,  3.92126485e-01],
       [ 2.01814509e-15,  2.28835628e-01,  3.92126485e-01]])

The calculator also support atom partitioned properties, including energies and charges.

>>> atoms.get_charges()  # in e
array([-0.66558478,  0.33279239,  0.33279239])
>>> atoms.get_potential_energies()  # in eV
array([-132.37935247,  -12.29407672,  -12.29407672])

ASE calculator#

class tblite.ase.TBLite(atoms: Atoms | None = None, **kwargs)[source]#

ASE calculator for using xTB Hamiltonians from the tblite library. Supported properties by this calculator are:

  • energy (free_energy): The total energy of the system (in eV)

  • energies: The atom-partitioned energy of the system (in eV)

  • forces: The derivative of the energy with respect to the atomic positions (in eV/Ang)

  • stress: The derivative of the energy with respect to cell deformations (PBC only, in eV/Ang³)

  • dipole: The dipole moment of the system (in Debye)

  • charges: Mulliken charges of the system (in e)

Supported keywords are

Keyword

Default

Description

method

“GFN2-xTB”

Underlying method for energy and forces

charge

None

Total charge of the system

multiplicity

None

Total multiplicity of the system

accuracy

1.0

Numerical accuracy of the calculation

electronic_temperature

300.0

Electronic temperatur in Kelvin

max_iterations

250

Iterations for self-consistent evaluation

initial_guess

“sad”

Initial guess for wavefunction (sad, eeq, or eeqbc)

mixer_damping

0.4

Damping parameter for self-consistent mixer

electric_field

None

Uniform electric field vector (in V/A)

spin_polarization

None

Spin polarization (scaling factor)

solvation

None

Solvation model to use (see below for details)

cache_api

True

Reuse generate API objects (recommended)

verbosity

1

Set verbosity of printout

xtb_config

None

Additional configuration for the API calculator (dict)

Solvation models are supported by passing a tuple to the solvation keyword. The first element of the tuple is the name of the solvation model, the following arguments are passed to initialize the solvation model.

Solvation model

Arguments

“alpb”

Solvent name (str), solution state (optional)

“gbsa”

Solvent name (str), solution state (optional)

“cpcm”

Epsilon (float)

“gbe”

Epsilon (float), Born kernel

“gb”

Epsilon (float), Born kernel

Example

An ASE calculator can be constructed by using the TBLite class provided by the tblite.ase module. For example to perform a single point calculation for a CO2 crystal use

>>> from tblite.ase import TBLite
>>> from ase.atoms import Atoms
>>> import numpy as np
>>> atoms = Atoms(
...     symbols="C4O8",
...     positions=np.array(
...         [
...             [0.9441259872, 0.9437851680, 0.9543505632],
...             [3.7179966528, 0.9556570368, 3.7316862240],
...             [3.7159517376, 3.7149292800, 0.9692330016],
...             [0.9529872864, 3.7220864832, 3.7296981120],
...             [1.6213905408, 1.6190616096, 1.6313879040],
...             [0.2656685664, 0.2694175776, 0.2776540416],
...             [4.3914553920, 1.6346256864, 3.0545920000],
...             [3.0440834880, 0.2764611744, 4.4080419264],
...             [4.3910577696, 3.0416409504, 0.2881058304],
...             [3.0399936576, 4.3879335936, 1.6497353376],
...             [0.2741322432, 4.4003734944, 3.0573754368],
...             [1.6312174944, 3.0434586528, 4.4023048032],
...         ]
...     ),
...     cell=np.array([5.68032, 5.68032, 5.68032]),
...     pbc=np.array([True, True, True]),
... )
>>> atoms.calc = TBLite(method="GFN1-xTB", verbosity=0)
>>> atoms.get_potential_energy()  # result in eV
-1257.0943962462964

The resulting calculator can be used like most ASE calculator, e.g. for optimizing geometries.

calculate(atoms: Atoms | None = None, properties: List[str] | None = None, system_changes: List[str] = ['positions', 'numbers', 'cell', 'pbc', 'initial_charges', 'initial_magmoms']) None[source]#

Perform actual calculation by calling the TBLite API

Example

>>> from ase.build import molecule
>>> from tblite.ase import TBLite
>>> calc = TBLite(method="GFN2-xTB")
>>> calc.calculate(molecule("H2O"))
>>> calc.get_potential_energy()
-137.96777625229421
>>> calc.calculate(molecule("CH4"))
>>> calc.get_potential_energy()
-113.60956621093894
Raises:
  • ase.calculators.calculator.InputError – on invalid input passed to the interface module

  • ase.calculators.calculator.CalculationFailed – in case of an RuntimeError in the library

default_parameters: dict[str, Any] = {'accuracy': 1.0, 'cache_api': True, 'charge': None, 'electric_field': None, 'electronic_temperature': 300.0, 'guess': 'sad', 'max_iterations': 250, 'method': 'GFN2-xTB', 'mixer_damping': 0.4, 'multiplicity': None, 'solvation': None, 'spin_polarization': None, 'verbosity': 1, 'xtb_config': None}#

Default parameters

implemented_properties: list[str] = ['energy', 'energies', 'forces', 'charges', 'dipole', 'stress']#

Properties calculator can handle (energy, forces, …)

reset() None[source]#

Clear all information from old calculation. This will only remove the cached API objects in case the cache_api is set to False.

set(**kwargs) dict[source]#

Set new parameters to TBLite. Will automatically reconstruct the underlying model in case critical parameters change.

Example

>>> from ase.build import molecule
>>> from tblite.ase import TBLite
>>> atoms = molecule("H2O")
>>> atoms.calc = TBLite(method="GFN2-xTB")
>>> atoms.get_potential_energy()
-137.96777625229421
>>> atoms.calc.set(method="GFN1-xTB")
{'method': 'GFN1-xTB'}
>>> atoms.get_potential_energy()
-156.9675057724589