Source code for ase.md.langevin

"""Langevin dynamics class."""

import numpy as np

from ase.md.md import MolecularDynamics
from ase.parallel import world, DummyMPI
from ase import units


[docs]class Langevin(MolecularDynamics): """Langevin (constant N, V, T) molecular dynamics.""" # Helps Asap doing the right thing. Increment when changing stuff: _lgv_version = 4 def __init__(self, atoms, timestep, temperature=None, friction=None, fixcm=True, *, temperature_K=None, trajectory=None, logfile=None, loginterval=1, communicator=world, rng=None, append_trajectory=False): """ Parameters: atoms: Atoms object The list of atoms. timestep: float The time step in ASE time units. temperature: float (deprecated) The desired temperature, in electron volt. temperature_K: float The desired temperature, in Kelvin. friction: float A friction coefficient, typically 1e-4 to 1e-2. fixcm: bool (optional) If True, the position and momentum of the center of mass is kept unperturbed. Default: True. rng: RNG object (optional) Random number generator, by default numpy.random. Must have a standard_normal method matching the signature of numpy.random.standard_normal. logfile: file object or str (optional) If *logfile* is a string, a file with that name will be opened. Use '-' for stdout. trajectory: Trajectory object or str (optional) Attach trajectory object. If *trajectory* is a string a Trajectory will be constructed. Use *None* (the default) for no trajectory. communicator: MPI communicator (optional) Communicator used to distribute random numbers to all tasks. Default: ase.parallel.world. Set to None to disable communication. append_trajectory: bool (optional) Defaults to False, which causes the trajectory file to be overwritten each time the dynamics is restarted from scratch. If True, the new structures are appended to the trajectory file instead. The temperature and friction are normally scalars, but in principle one quantity per atom could be specified by giving an array. RATTLE constraints can be used with these propagators, see: E. V.-Eijnden, and G. Ciccotti, Chem. Phys. Lett. 429, 310 (2006) The propagator is Equation 23 (Eq. 39 if RATTLE constraints are used) of the above reference. That reference also contains another propagator in Eq. 21/34; but that propagator is not quasi-symplectic and gives a systematic offset in the temperature at large time steps. """ if friction is None: raise TypeError("Missing 'friction' argument.") self.fr = friction self.temp = units.kB * self._process_temperature(temperature, temperature_K, 'eV') self.fix_com = fixcm if communicator is None: communicator = DummyMPI() self.communicator = communicator if rng is None: self.rng = np.random else: self.rng = rng MolecularDynamics.__init__(self, atoms, timestep, trajectory, logfile, loginterval, append_trajectory=append_trajectory) self.updatevars() def todict(self): d = MolecularDynamics.todict(self) d.update({'temperature_K': self.temp / units.kB, 'friction': self.fr, 'fixcm': self.fix_com}) return d def set_temperature(self, temperature=None, temperature_K=None): self.temp = units.kB * self._process_temperature(temperature, temperature_K, 'eV') self.updatevars() def set_friction(self, friction): self.fr = friction self.updatevars() def set_timestep(self, timestep): self.dt = timestep self.updatevars() def updatevars(self): dt = self.dt T = self.temp fr = self.fr masses = self.masses sigma = np.sqrt(2 * T * fr / masses) self.c1 = dt / 2. - dt * dt * fr / 8. self.c2 = dt * fr / 2 - dt * dt * fr * fr / 8. self.c3 = np.sqrt(dt) * sigma / 2. - dt**1.5 * fr * sigma / 8. self.c5 = dt**1.5 * sigma / (2 * np.sqrt(3)) self.c4 = fr / 2. * self.c5 def step(self, forces=None): atoms = self.atoms natoms = len(atoms) if forces is None: forces = atoms.get_forces(md=True) # This velocity as well as xi, eta and a few other variables are stored # as attributes, so Asap can do its magic when atoms migrate between # processors. self.v = atoms.get_velocities() self.xi = self.rng.standard_normal(size=(natoms, 3)) self.eta = self.rng.standard_normal(size=(natoms, 3)) # When holonomic constraints for rigid linear triatomic molecules are # present, ask the constraints to redistribute xi and eta within each # triple defined in the constraints. This is needed to achieve the # correct target temperature. for constraint in self.atoms.constraints: if hasattr(constraint, 'redistribute_forces_md'): constraint.redistribute_forces_md(atoms, self.xi, rand=True) constraint.redistribute_forces_md(atoms, self.eta, rand=True) self.communicator.broadcast(self.xi, 0) self.communicator.broadcast(self.eta, 0) # First halfstep in the velocity. self.v += (self.c1 * forces / self.masses - self.c2 * self.v + self.c3 * self.xi - self.c4 * self.eta) # Full step in positions x = atoms.get_positions() if self.fix_com: old_com = atoms.get_center_of_mass() # Step: x^n -> x^(n+1) - this applies constraints if any. atoms.set_positions(x + self.dt * self.v + self.c5 * self.eta) if self.fix_com: atoms.set_center_of_mass(old_com) # recalc velocities after RATTLE constraints are applied self.v = (self.atoms.get_positions() - x - self.c5 * self.eta) / self.dt forces = atoms.get_forces(md=True) # Update the velocities self.v += (self.c1 * forces / self.masses - self.c2 * self.v + self.c3 * self.xi - self.c4 * self.eta) if self.fix_com: # subtract center of mass vel self.v -= self._get_com_velocity(self.v) # Second part of RATTLE taken care of here atoms.set_momenta(self.v * self.masses) return forces