1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108
|
.. module:: ase.calculators.qmmm
QMMM
====
There are three QM/MM calculators native to ASE:
========================= ===================
Explicit Interaction QMMM :class:`EIQMMM`
Simple, subtrative QMMM :class:`SimpleQMMM`
Force-based qmmm :class:`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
:doi:`in this paper <10.1021/acs.jctc.7b00621>`, which should be
cited if the method is used.
Other ASE-calculators that currently support EIQMMM:
1. :mod:`DFTB+ <ase.calculators.dftb>`
2. :mod:`CRYSTAL <ase.calculators.crystal>`
3. :mod:`TURBOMOLE <ase.calculators.turbomole>`
4. :mod:`ORCA <ase.calculators.orca>`
.. _GPAW: https://gpaw.readthedocs.io
.. seealso::
The :ref:`qmmm` tutorial, on how to use the Explicit Interaction QMMM
calculator
.. autoclass:: EIQMMM
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:
.. autoclass:: LJInteractions
Or, for couplings requiring more generality, you should use:
.. autoclass:: LJInteractionsGeneral
You can control how the QM part is embedded in the MM part by supplying your
own embedding object when you construct the :class:`EIQMMM` instance. The
Embedding object will be specific to the QM calculator you want to use. The
default is this one:
.. autoclass:: 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.
.. autoclass:: SimpleQMMM
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
:doi:`N. Bernstein, J. R. Kermode, and G. Csányi, Rep. Prog. Phys. 72, 026501 (2009) <10.1088/0034-4885/72/2/026501>`
for a review of force-based QM/MM approaches, which should be cited if this method is used,
and :doi:`T. D. Swinburne and J. R. Kermode, Phys. Rev. B 96, 144102 (2017) <10.1103/PhysRevB.96.144102>`
for an application which used this implementation.
.. autoclass:: ForceQMMM
Basic usage is as follows::
from ase.calculators.qmmm import ForceQMMM
atoms = ...
qm_mask = ...
atoms.calc = ForceQMMM(qm_mask,
QMCalculator(...),
MMCalculator(...),
buffer_width=...)
See :git:`ase/test/forcefields/test_forceqmmm.py` test-case for a complete example.
|