QMMM

There are three QM/MM calculators native to ASE:

Explicit Interaction QMMM

EIQMMM

Simple, subtrative QMMM

SimpleQMMM

Force-based qmmm

ForceQMMM

Explicit Interaction QMMM

In Explicit Interaction QMMM, the QM and MM regions are explicitly coupled with an electrostatic interaction term. This requires that the electrostatic potential from the classical charges of the MM subsystem is fed into the QM calculator. This is built into GPAW. More info in this paper, which should be cited if the method is used.

Other ASE-calculators that currently support EIQMMM:

  1. DFTB+

  2. CRYSTAL

  3. TURBOMOLE

  4. ORCA

See also

The ASE for QM/MM Simulations tutorial, on how to use the Explicit Interaction QMMM calculator

class ase.calculators.qmmm.EIQMMM(selection, qmcalc, mmcalc, interaction, vacuum=None, embedding=None, output=None)[source]

Explicit interaction QMMM calculator.

EIQMMM object.

The energy is calculated as:

        _          _         _    _
E = E  (R  ) + E  (R  ) + E (R  , R  )
     QM  QM     MM  MM     I  QM   MM

parameters:

selection: list of int, slice object or list of bool

Selection out of all the atoms that belong to the QM part.

qmcalc: Calculator object

QM-calculator.

mmcalc: Calculator object

MM-calculator.

interaction: Interaction object

Interaction between QM and MM regions.

vacuum: float or None

Amount of vacuum to add around QM atoms. Use None if QM calculator doesn’t need a box.

embedding: Embedding object or None

Specialized embedding object. Use None in order to use the default one.

output: None, ‘-’, str or file-descriptor.

File for logging information - default is no logging (None).

Here, you need to specify the interaction:

from ase.calculators.qmmm import EIQMMM, LJInteraction
from ase.calculators.tip3p import epsilon0, sigma0
lj = LJInteraction({'OO': (epsilon0, sigma0)})
atoms.calc = EIQMMM([0, 1, 2],
                    QMCalculator(...),
                    MMCalculator(...),
                    interaction=lj)

For Lennard-Jones type of interactions you can use:

class ase.calculators.qmmm.LJInteractions(parameters)[source]

Lennard-Jones type explicit interaction.

parameters: dict

Mapping from pair of atoms to tuple containing epsilon and sigma for that pair.

Example

lj = LJInteractions({(‘O’, ‘O’): (eps, sigma)})

Or, for couplings requiring more generality, you should use:

class ase.calculators.qmmm.LJInteractionsGeneral(sigmaqm, epsilonqm, sigmamm, epsilonmm, qm_molecule_size, mm_molecule_size=3, rc=inf, width=1.0)[source]

General Lennard-Jones type explicit interaction.

sigmaqm: array

Array of sigma-parameters which should have the length of the QM subsystem

epsilonqm: array

As sigmaqm, but for epsilon-paramaters

sigmamm: Either array (A) or tuple (B)
A (no counterions):

Array of sigma-parameters with the length of the smallests repeating atoms-group (i.e. molecule) of the MM subsystem

B (counterions):

Tuple: (arr1, arr2), where arr1 is an array of sigmas with the length of counterions in the MM subsystem, and arr2 is the array from A.

epsilonmm: array or tuple

As sigmamm but for epsilon-parameters.

qm_molecule_size: int

number of atoms of the smallest repeating atoms-group (i.e. molecule) in the QM subsystem (often just the number of atoms in the QM subsystem)

mm_molecule_size: int

as qm_molecule_size but for the MM subsystem. Will be overwritten if counterions are present in the MM subsystem (via the CombineMM calculator)

You can control how the QM part is embedded in the MM part by supplying your own embedding object when you construct the EIQMMM instance. The Embedding object will be specific to the QM calculator you want to use. The default is this one:

class ase.calculators.qmmm.Embedding(molecule_size=3, **parameters)[source]

Point-charge embedding.

Simple, subtractive QMMM calculations

This QM/MM calculator is similar to the original ONIOM model, doing simple, subtractive QM/MM between any two calculators.

class ase.calculators.qmmm.SimpleQMMM(selection, qmcalc, mmcalc1, mmcalc2, vacuum=None)[source]

Simple QMMM calculator.

SimpleQMMM object.

The energy is calculated as:

        _          _          _
E = E  (R  ) - E  (R  ) + E  (R   )
     QM  QM     MM  QM     MM  all

parameters:

selection: list of int, slice object or list of bool

Selection out of all the atoms that belong to the QM part.

qmcalc: Calculator object

QM-calculator.

mmcalc1: Calculator object

MM-calculator used for QM region.

mmcalc2: Calculator object

MM-calculator used for everything.

vacuum: float or None

Amount of vacuum to add around QM atoms. Use None if QM calculator doesn’t need a box.

This type of QMMM can combine any pair of ASE calculators:

from ase.calculators.qmmm import SimpleQMMM
atoms = ...
atoms.calc = SimpleQMMM([0, 1, 2],
                        QMCalculator(...),
                        MMCalculator(...))

where [0, 1, 2] would select the first three atoms for the QM-part.

Force-based QM/MM

This QM/MM calculator mixes forces from any pair of ASE calculators. A finite buffer is added around the core QM region to ensure accurate forces; careful testing of the required buffer size is required. See N. Bernstein, J. R. Kermode, and G. Csányi, Rep. Prog. Phys. 72, 026501 (2009) for a review of force-based QM/MM approaches, which should be cited if this method is used, and T. D. Swinburne and J. R. Kermode, Phys. Rev. B 96, 144102 (2017) for an application which used this implementation.

class ase.calculators.qmmm.ForceQMMM(atoms, qm_selection_mask, qm_calc, mm_calc, buffer_width, vacuum=5.0, zero_mean=True, qm_cell_round_off=3, qm_radius=None)[source]

Force-based QM/MM calculator

QM forces are computed using a buffer region and then mixed abruptly with MM forces:

F^i_QMMM = { F^i_QM if i in QM region

{ F^i_MM otherwise

cf. N. Bernstein, J. R. Kermode, and G. Csanyi, Rep. Prog. Phys. 72, 026501 (2009) and T. D. Swinburne and J. R. Kermode, Phys. Rev. B 96, 144102 (2017).

ForceQMMM calculator

Parameters:

qm_selection_mask: list of ints, slice object or bool list/array

Selection out of atoms that belong to the QM region.

qm_calc: Calculator object

QM-calculator.

mm_calc: Calculator object

MM-calculator (should be scaled, see RescaledCalculator) Can use \(ForceConstantCalculator\) based on QM force constants, if available.

vacuum: float or None

Amount of vacuum to add around QM atoms.

zero_mean: bool

If True, add a correction to zero the mean force in each direction

qm_cell_round_off: float

Tolerance value in Angstrom to round the qm cluster cell

qm_radius: 3x1 array of floats qm_radius for [x, y, z]

3d qm radius for calculation of qm cluster cell. default is None and the radius is estimated from maximum distance between the atoms in qm region.

Basic usage is as follows:

from ase.calculators.qmmm import ForceQMMM
atoms = ...
qm_mask = ...
atoms.calc = ForceQMMM(qm_mask,
                       QMCalculator(...),
                       MMCalculator(...),
                       buffer_width=...)

See ase/test/forcefields/test_forceqmmm.py test-case for a complete example.