Calculators

For ASE, a calculator is a black box that can take atomic numbers and atomic positions from an Atoms object and calculate the energy and forces and sometimes also stresses.

In order to calculate forces and energies, you need to attach a calculator object to your atoms object:

>>> atoms = read('molecule.xyz')
>>> e = atoms.get_potential_energy()  
Traceback (most recent call last):
  File "<stdin>", line 1, in <module>
  File "/home/jjmo/ase/atoms/ase.py", line 399, in get_potential_energy
    raise RuntimeError('Atoms object has no calculator.')
RuntimeError: Atoms object has no calculator.
>>> from ase.calculators.abinit import Abinit
>>> calc = Abinit(...)
>>> atoms.calc = calc
>>> e = atoms.get_potential_energy()
>>> print(e)
-42.0

Here we attached an instance of the ase.calculators.abinit class and then we asked for the energy.

Supported calculators

The calculators can be divided in four groups:

1) Asap, BigDFT, DFTK, GPAW, and Hotbit have their own native ASE interfaces.

  1. ABINIT, AMBER, CP2K, CASTEP, deMon2k, DFTB+, ELK, EXCITING, FHI-aims, FLEUR, GAUSSIAN, Gromacs, LAMMPS, MOPAC, NWChem, Octopus, ONETEP, psi4, Q-Chem, Quantum ESPRESSO, SIESTA, TURBOMOLE and VASP, have Python wrappers in the ASE package, but the actual FORTRAN/C/C++ codes are not part of ASE.

  2. Pure python implementations included in the ASE package: EMT, EAM, Lennard-Jones and Morse.

  3. Calculators that wrap others, included in the ASE package: ase.calculators.checkpoint.CheckpointCalculator, the ase.calculators.loggingcalc.LoggingCalculator, the ase.calculators.socketio.SocketIOCalculator, the Grimme-D3 potential, and the qmmm calculators EIQMMM, and SimpleQMMM.

name

description

Asap

Highly efficient EMT code

BigDFT

Wavelet based code for DFT

DFTK

Plane-wave code for DFT and related models

GPAW

Real-space/plane-wave/LCAO PAW code

Hotbit

DFT based tight binding

abinit

Plane-wave pseudopotential code

amber

Classical molecular dynamics code

castep

Plane-wave pseudopotential code

cp2k

DFT and classical potentials

demon

Gaussian based DFT code

demonnano

DFT based tight binding code

dftb

DFT based tight binding

dmol

Atomic orbital DFT code

eam

Embedded Atom Method

elk

Full Potential LAPW code

espresso

Plane-wave pseudopotential code

exciting

Full Potential LAPW code

aims

Numeric atomic orbital, full potential code

fleur

Full Potential LAPW code

gamess_us

Gaussian based electronic structure code

gaussian

Gaussian based electronic structure code

gromacs

Classical molecular dynamics code

gulp

Interatomic potential code

kim

Classical MD with standardized models

lammps

Classical molecular dynamics code

mixing

Combination of multiple calculators

mopac

Semiempirical molecular orbital code

nwchem

Gaussian based electronic structure code

octopus

Real-space pseudopotential code

onetep

Linear-scaling pseudopotential code

openmx

LCAO pseudopotential code

orca

Gaussian based electronic structure code

psi4

Gaussian based electronic structure code

qchem

Gaussian based electronic structure code

siesta

LCAO pseudopotential code

turbomole

Fast atom orbital code

vasp

Plane-wave PAW code

emt

Effective Medium Theory calculator

lj

Lennard-Jones potential

morse

Morse potential

checkpoint

Checkpoint calculator

socketio

Socket-based interface to calculators

loggingcalc

Logging calculator

dftd3

DFT-D3 dispersion correction calculator

EIQMMM

Explicit Interaction QM/MM

SimpleQMMM

Subtractive (ONIOM style) QM/MM

Note

A Fortran implemetation of the Grimme-D3 potential, that can be used as an add-on to any ASE calculator, can be found here: https://gitlab.com/ehermes/ased3/tree/master.

The calculators included in ASE are used like this:

>>> from ase.calculators.abc import ABC
>>> calc = ABC(...)

where abc is the module name and ABC is the class name.

Calculator keywords

Example for a hypothetical ABC calculator:

ABC(restart=None, ignore_bad_restart_file=False, label=None,
atoms=None, parameters=None, command='abc > PREFIX.abc',
xc=None, kpts=[1, 1, 1], smearing=None,
charge=0.0, nbands=None, **kwargs)

Create ABC calculator

restart: str

Prefix for restart file. May contain a directory. Default is None: don’t restart.

ignore_bad_restart_file: bool

Ignore broken or missing restart file. By default, it is an error if the restart file is missing or broken.

label: str

Name used for all files. May contain a directory.

atoms: Atoms object

Optional Atoms object to which the calculator will be attached. When restarting, atoms will get its positions and unit-cell updated from file.

command: str

Command used to start calculation. This will override any value in an ASE_ABC_COMMAND environment variable.

parameters: str

Read parameters from file.

xc: str

XC-functional ('LDA', 'PBE', …).

kpts:

Brillouin zone sampling:

  • (1,1,1): Gamma-point

  • (n1,n2,n3): Monkhorst-Pack grid

  • (n1,n2,n3,'gamma'): Shifted Monkhorst-Pack grid that includes \(\Gamma\)

  • [(k11,k12,k13),(k21,k22,k23),...]: Explicit list in units of the reciprocal lattice vectors

  • kpts=3.5: \(\vec k\)-point density as in 3.5 \(\vec k\)-points per Å\(^{-1}\).

smearing: tuple

The smearing of occupation numbers. Must be a tuple:

  • ('Fermi-Dirac', width)

  • ('Gaussian', width)

  • ('Methfessel-Paxton', width, n), where \(n\) is the order (\(n=0\) is the same as 'Gaussian')

Lower-case names are also allowed. The width parameter is given in eV units.

charge: float

Charge of the system in units of \(|e|\) (charge=1 means one electron has been removed). Default is charge=0.

nbands: int

Number of bands. Each band can be occupied by two electrons.

Not all of the above arguments make sense for all of ASE’s calculators. As an example, Gromacs will not accept DFT related keywords such as xc and smearing. In addition to the keywords mentioned above, each calculator may have native keywords that are specific to only that calculator.

Keyword arguments can also be set or changed at a later stage using the set() method:

ase.calculators.set(key1=value1, key2=value2, ...)

Calculator interface

All calculators must have the following interface:

class ase.calculators.interface.Calculator[source]

ASE calculator.

A calculator should store a copy of the atoms object used for the last calculation. When one of the get_potential_energy, get_forces, or get_stress methods is called, the calculator should check if anything has changed since the last calculation and only do the calculation if it’s really needed. Two sets of atoms are considered identical if they have the same positions, atomic numbers, unit cell and periodic boundary conditions.

calculation_required(atoms, quantities)[source]

Check if a calculation is required.

Check if the quantities in the quantities list have already been calculated for the atomic configuration atoms. The quantities can be one or more of: ‘energy’, ‘forces’, ‘stress’, ‘charges’ and ‘magmoms’.

This method is used to check if a quantity is available without further calculations. For this reason, calculators should react to unknown/unsupported quantities by returning True, indicating that the quantity is not available.

get_forces(atoms)[source]

Return the forces.

get_potential_energy(atoms=None, force_consistent=False)[source]

Return total energy.

Both the energy extrapolated to zero Kelvin and the energy consistent with the forces (the free energy) can be returned.

get_stress(atoms)[source]

Return the stress.

Electronic structure calculators

These calculators have wave functions, electron densities, eigenvalues and many other quantities. Therefore, it makes sense to have a set of standard methods for accessing those quantities:

class ase.calculators.interface.DFTCalculator[source]

Class for demonstrating the ASE interface to DFT-calculators.

get_bz_k_points()[source]

Return all the k-points in the 1. Brillouin zone.

The coordinates are relative to reciprocal latice vectors.

get_effective_potential(spin=0, pad=True)[source]

Return pseudo-effective-potential array.

get_eigenvalues(kpt=0, spin=0)[source]

Return eigenvalue array.

get_fermi_level()[source]

Return the Fermi level.

get_ibz_k_points()[source]

Return k-points in the irreducible part of the Brillouin zone.

The coordinates are relative to reciprocal latice vectors.

get_k_point_weights()[source]

Weights of the k-points.

The sum of all weights is one.

get_magnetic_moment(atoms=None)[source]

Return the total magnetic moment.

get_number_of_bands()[source]

Return the number of bands.

get_number_of_grid_points()[source]

Return the shape of arrays.

get_number_of_spins()[source]

Return the number of spins in the calculation.

Spin-paired calculations: 1, spin-polarized calculation: 2.

get_occupation_numbers(kpt=0, spin=0)[source]

Return occupation number array.

get_pseudo_density(spin=None, pad=True)[source]

Return pseudo-density array.

If spin is not given, then the total density is returned. Otherwise, the spin up or down density is returned (spin=0 or 1).

get_pseudo_wave_function(band=0, kpt=0, spin=0, broadcast=True, pad=True)[source]

Return pseudo-wave-function array.

get_spin_polarized()[source]

Is it a spin-polarized calculation?

get_wannier_localization_matrix(nbands, dirG, kpoint, nextkpoint, G_I, spin)[source]

Calculate integrals for maximally localized Wannier functions.

get_xc_functional()[source]

Return the XC-functional identifier.

‘LDA’, ‘PBE’, …

initial_wannier(initialwannier, kpointgrid, fixedstates, edf, spin, nbands)[source]

Initial guess for the shape of wannier functions.

Use initial guess for wannier orbitals to determine rotation matrices U and C.