Dimer method
The dimer method is a means of finding a saddle point on a potential energy surface starting from a single point (as opposed to a NEB calculation, which requires an initial and final state). You can read about this method here:
‘A dimer method for finding saddle points on high dimensional potential surfaces using only first derivatives’, G. Henkelman and H. Jonsson, J. Chem. Phys. 111, 7010, 1999.
An example is shown below.
import pytest
from ase.build import add_adsorbate, fcc100
from ase.calculators.emt import EMT
from ase.constraints import FixAtoms
from ase.mep import DimerControl, MinModeAtoms, MinModeTranslate
@pytest.mark.optimize()
def test_dimer_method(testdir):
# Set up a small "slab" with an adatoms
atoms = fcc100('Pt', size=(2, 2, 1), vacuum=10.0)
add_adsorbate(atoms, 'Pt', 1.611, 'hollow')
# Freeze the "slab"
mask = [atom.tag > 0 for atom in atoms]
atoms.set_constraint(FixAtoms(mask=mask))
# Calculate using EMT
atoms.calc = EMT()
atoms.get_potential_energy()
# Set up the dimer
with DimerControl(initial_eigenmode_method='displacement',
displacement_method='vector', logfile=None,
mask=[0, 0, 0, 0, 1]) as d_control:
d_atoms = MinModeAtoms(atoms, d_control)
# Displace the atoms
displacement_vector = [[0.0] * 3] * 5
displacement_vector[-1][1] = -0.1
d_atoms.displace(displacement_vector=displacement_vector)
# Converge to a saddle point
with MinModeTranslate(d_atoms, trajectory='dimer_method.traj',
logfile=None) as dim_rlx:
dim_rlx.run(fmax=0.001)
The module contains several classes.
- class ase.mep.dimer.DimerControl(logfile='-', eigenmode_logfile=None, comm=<ase.parallel.MPI object>, **kwargs)[source]
A class that takes care of the parameters needed for a Dimer search.
Parameters:
- eigenmode_method: str
The name of the eigenmode search method.
- f_rot_min: float
Size of the rotational force under which no rotation will be performed.
- f_rot_max: float
Size of the rotational force under which only one rotation will be performed.
- max_num_rot: int
Maximum number of rotations per optimizer step.
- trial_angle: float
Trial angle for the finite difference estimate of the rotational angle in radians.
- trial_trans_step: float
Trial step size for the MinModeTranslate optimizer.
- maximum_translation: float
Maximum step size and forced step size when the curvature is still positive for the MinModeTranslate optimizer.
- cg_translation: bool
Conjugate Gradient for the MinModeTranslate optimizer.
- use_central_forces: bool
Only calculate the forces at one end of the dimer and extrapolate the forces to the other.
- dimer_separation: float
Separation of the dimer’s images.
- initial_eigenmode_method: str
How to construct the initial eigenmode of the dimer. If an eigenmode is given when creating the MinModeAtoms object, this will be ignored. Possible choices are: ‘gauss’ and ‘displacement’
- extrapolate_forces: bool
When more than one rotation is performed, an extrapolation scheme can be used to reduce the number of force evaluations.
- displacement_method: str
How to displace the atoms. Possible choices are ‘gauss’ and ‘vector’.
- gauss_std: float
The standard deviation of the gauss curve used when doing random displacement.
- order: int
How many lowest eigenmodes will be inverted.
- mask: list of bool
Which atoms will be moved during displacement.
- displacement_center: int or [float, float, float]
The center of displacement, nearby atoms will be displaced.
- displacement_radius: float
When choosing which atoms to displace with the displacement_center keyword, this decides how many nearby atoms to displace.
- number_of_displacement_atoms: int
The amount of atoms near displacement_center to displace.
- class ase.mep.dimer.MinModeAtoms(atoms, control=None, eigenmodes=None, random_seed=None, comm=<ase.parallel.MPI object>, **kwargs)[source]
Wrapper for Atoms with information related to minimum mode searching.
Contains an Atoms object and pipes all unknown function calls to that object. Other information that is stored in this object are the estimate for the lowest eigenvalue, curvature, and its corresponding eigenmode, eigenmode. Furthermore, the original configuration of the Atoms object is stored for use in multiple minimum mode searches. The forces on the system are modified by inverting the component along the eigenmode estimate. This eventually brings the system to a saddle point.
Parameters:
- atomsAtoms object
A regular Atoms object
- controlMinModeControl object
Contains the parameters necessary for the eigenmode search. If no control object is supplied a default DimerControl will be created and used.
- mask: list of bool
Determines which atoms will be moved when calling displace()
- random_seed: int
The seed used for the random number generator. Defaults to modified version the current time.
- class ase.mep.dimer.MinModeTranslate(dimeratoms, logfile='-', trajectory=None)[source]
An Optimizer specifically tailored to minimum mode following.
- Parameters:
atoms (
Atoms
) – The Atoms object to relax.restart (str) – Filename for restart file. Default value is None.
logfile (file object or str) – If logfile is a string, a file with that name will be opened. Use ‘-’ for stdout.
trajectory (Trajectory object or str) – Attach trajectory object. If trajectory is a string a Trajectory will be constructed. Use None for no trajectory.
append_trajectory (bool) – Appended to the trajectory file instead of overwriting it.
kwargs (dict, optional) – Extra arguments passed to
Dynamics
.
- class ase.mep.dimer.DimerEigenmodeSearch(dimeratoms, control=None, eigenmode=None, basis=None, **kwargs)[source]
An implementation of the Dimer’s minimum eigenvalue mode search.
This class implements the rotational part of the dimer saddle point searching method.
Parameters:
- atoms: MinModeAtoms object
MinModeAtoms is an extension to the Atoms object, which includes information about the lowest eigenvalue mode.
- control: DimerControl object
Contains the parameters necessary for the eigenmode search. If no control object is supplied a default DimerControl will be created and used.
- basis: list of xyz-values
Eigenmode. Must be an ndarray of shape (n, 3). It is possible to constrain the eigenmodes to be orthogonal to this given eigenmode.
Notes:
The code is inspired, with permission, by code written by the Henkelman group, which can be found at http://theory.cm.utexas.edu/vtsttools/code/
References:
Henkelman and Jonsson, JCP 111, 7010 (1999)
Olsen, Kroes, Henkelman, Arnaldsson, and Jonsson, JCP 121, 9776 (2004).
Heyden, Bell, and Keil, JCP 123, 224101 (2005).
Kastner and Sherwood, JCP 128, 014106 (2008).