integrators¶
-
class
atomsmm.integrators.
GlobalThermostatIntegrator
(stepSize, nveIntegrator, thermostat=None)[source]¶ This class extends OpenMM’s CustomIntegrator class in order to facilitate the construction of NVT integrators which include a global thermostat, that is, one that acts equally and simultaneously on all degrees of freedom of the system. In this case, a complete NVT step is split as:
\[e^{\delta t \, iL_\mathrm{NVT}} = e^{\frac{1}{2} \delta t \, iL_\mathrm{T}} e^{\delta t \, iL_\mathrm{NVE}} e^{\frac{1}{2} \delta t \, iL_\mathrm{T}}\]The propagator \(e^{\delta t \, iL_\mathrm{NVE}}\) is a Hamiltonian
corresponds to a Hamiltonian \(iL_\mathrm{T}\)
- Parameters
stepSize (unit.Quantity) – The step size with which to integrate the system (in time unit).
nveIntegrator (
HamiltonianPropagator
) – The Hamiltonian propagator.thermostat (
ThermostatPropagator
, optional, default=None) – The thermostat propagator.randomSeed (int, optional, default=None) – A seed for random numbers.
-
class
atomsmm.integrators.
MultipleTimeScaleIntegrator
(stepSize, loops, move=None, boost=None, bath=None, **kwargs)[source]¶ This class implements a Multiple Time-Scale (MTS) integrator using the RESPA method.
- Parameters
stepSize (unit.Quantity) – The largest time step for numerically integrating the system of equations.
loops (list(int)) – A list of N integers. Assuming that force group 0 is composed of the fastest forces, while group N-1 is composed of the slowest ones, loops[k] determines how many steps involving forces of group k are internally executed for every step involving those of group k+1.
move (
Propagator
, optional, default = None) – A move propagator.boost (
Propagator
, optional, default = None) – A boost propagator.bath (
Propagator
, optional, default = None) – A bath propagator.
- Keyword Arguments
scheme (str, optional, default = middle) – The splitting scheme used to solve the equations of motion. Available options are middle, xi-respa, xo-respa, side, and blitz. If it is middle (default), then the bath propagator will be inserted between half-step coordinate moves during the fastest-force loops. If it is xi-respa, xo-respa, or side, then the bath propagator will be integrated in both extremities of each loop concerning one of the N time scales, with xi-respa referring to the time scale of fastest forces (force group 0), xo-respa referring to the time scale of the slowest forces (force group N-1), and side requiring the user to select the time scale in which to locate the bath propagator via keyword argument location (see below). If it is blitz, then the force-related propagators will be fully integrated at the outset of each loop in all time scales and the bath propagator will be integrated between half-step coordinate moves during the fastest-force loops.
location (int, optional, default = None) – The index of the force group (from 0 to N-1) that defines the time scale in which the bath propagator will be located. This is only meaningful if keyword scheme is set to side (see above).
nsy (int, optional, default = 1) – The number of Suzuki-Yoshida terms to factorize the bath propagator. Valid options are 1, 3, 7, and 15.
nres (int, optional, default = 1) – The number of RESPA-like subdivisions to factorize the bath propagator.
Warning
The xo-respa and xi-respa schemes implemented here are slightly different from the ones described in the paper by Leimkuhler, Margul, and Tuckerman [4].
-
class
atomsmm.integrators.
NHL_R_Integrator
(stepSize, loops, temperature, timeScale, frictionConstant, **kwargs)[source]¶ This class is an implementation of the massive Nosé-Hoover-Langevin (RESPA) integrator. The method consists in solving the following equations for each degree of freedom (DOF) in the system:
\[\begin{split}& \frac{dx}{dt} = v \\ & \frac{dv}{dt} = \frac{f}{m} - v_2 v \\ & dv_2 = \frac{m v^2 - kT}{Q_2}dt - \gamma v_2 dt + \sqrt{\frac{2 \gamma kT}{Q_2}} dW\end{split}\]The equations are integrated by a reversible, multiple timescale numerical scheme.
- Parameters
stepSize (unit.Quantity) – The largest time step for numerically integrating the system of equations.
loops (list(int)) – See description in
MultipleTimeScaleIntegrator
.temperature (unit.Quantity) – The temperature to which the configurational sampling should correspond.
timeScale (unit.Quantity) – A time scale \(\tau\) from which the inertial parameters are computed as \(Q_2 = kT\tau^2\).
frictionConstant (unit.Quantity) – The friction constant \(\gamma\) present in the stochastic equation of motion for per-DOF thermostat variable \(v_2\).
**kwargs (keyword arguments) – The same keyword arguments of class
MultipleTimeScaleIntegrator
apply here.
-
class
atomsmm.integrators.
Langevin_R_Integrator
(stepSize, loops, temperature, frictionConstant, **kwargs)[source]¶ This class is an implementation of the multiple time scale Langevin (RESPA) integrator. The method consists in solving the following equations for each degree of freedom (DOF) in the system:
\[\begin{split}& \frac{dx}{dt} = v \\ & \frac{dv}{dt} = \frac{f}{m} - \gamma v dt + \sqrt{\frac{2 \gamma kT}{m}} dW\end{split}\]The equations are integrated by a reversible, multiple timescale numerical scheme.
- Parameters
stepSize (unit.Quantity) – The largest time step for numerically integrating the system of equations.
loops (list(int)) – See description in
MultipleTimeScaleIntegrator
.temperature (unit.Quantity) – The temperature to which the configurational sampling should correspond.
frictionConstant (unit.Quantity) – The friction constant \(\gamma\) present in the stochastic equation of motion for per-DOF thermostat variable \(v_2\).
**kwargs (keyword arguments) – The same keyword arguments of class
MultipleTimeScaleIntegrator
apply here.
-
class
atomsmm.integrators.
SIN_R_Integrator
(stepSize, loops, temperature, timeScale, frictionConstant, **kwargs)[source]¶ This class is an implementation of the Stochastic-Iso-NH-RESPA or SIN(R) method of Leimkuhler, Margul, and Tuckerman [4]. The method consists in solving the following equations for each degree of freedom (DOF) in the system:
\[\begin{split}& \frac{dx}{dt} = v \\ & \frac{dv}{dt} = \frac{f}{m} - \lambda v \\ & \frac{dv_1}{dt} = - \lambda v_1 - v_2 v_1 \\ & dv_2 = \frac{Q_1 v_1^2 - kT}{Q_2}dt - \gamma v_2 dt + \sqrt{\frac{2 \gamma kT}{Q_2}} dW\end{split}\]where:
\[\lambda = \frac{f v - \frac{1}{2} Q_1 v_2 v_1^2}{m v^2 + \frac{1}{2} Q_1 v_1^2}.\]A consequence of these equations is that
\[m v^2 + \frac{1}{2} Q_1 v_1^2 = kT.\]The equations are integrated by a reversible, multiple timescale numerical scheme.
- Parameters
stepSize (unit.Quantity) – The largest time step for numerically integrating the system of equations.
loops (list(int)) – See description in
MultipleTimeScaleIntegrator
.temperature (unit.Quantity) – The temperature to which the configurational sampling should correspond.
timeScale (unit.Quantity) – A time scale \(\tau\) from which the inertial parameters are computed as \(Q_1 = Q_2 = kT\tau^2\).
frictionConstant (unit.Quantity) – The friction constant \(\gamma\) present in the stochastic equation of motion for per-DOF thermostat variable \(v_2\).
**kwargs (keyword arguments) – The same keyword arguments of class
MultipleTimeScaleIntegrator
apply here.
-
class
atomsmm.integrators.
NewMethodIntegrator
(stepSize, loops, temperature, timeScale, frictionConstant, **kwargs)[source]¶ This class is an implementation of the Stochastic-Iso-NH-RESPA or SIN(R) method of Leimkuhler, Margul, and Tuckerman [4]. The method consists in solving the following equations for each degree of freedom (DOF) in the system:
\[\begin{split}& \frac{dx}{dt} = v \\ & \frac{dv}{dt} = \frac{f}{m} - \lambda v \\ & \frac{dv_1}{dt} = - \lambda v_1 - v_2 v_1 \\ & dv_2 = \frac{Q_1 v_1^2 - kT}{Q_2}dt - \gamma v_2 dt + \sqrt{\frac{2 \gamma kT}{Q_2}} dW\end{split}\]where:
\[\lambda = \frac{f v - \frac{1}{2} Q_1 v_2 v_1^2}{m v^2 + \frac{1}{2} Q_1 v_1^2}.\]A consequence of these equations is that
\[m v^2 + \frac{1}{2} Q_1 v_1^2 = kT.\]The equations are integrated by a reversible, multiple timescale numerical scheme.
- Parameters
stepSize (unit.Quantity) – The largest time step for numerically integrating the system of equations.
loops (list(int)) – See description in
MultipleTimeScaleIntegrator
.temperature (unit.Quantity) – The temperature to which the configurational sampling should correspond.
timeScale (unit.Quantity) – A time scale \(\tau\) from which the inertial parameters are computed as \(Q_1 = Q_2 = kT\tau^2\).
frictionConstant (unit.Quantity) – The friction constant \(\gamma\) present in the stochastic equation of motion for per-DOF thermostat variable \(v_2\).
**kwargs (keyword arguments) – The same keyword arguments of class
MultipleTimeScaleIntegrator
apply here.
-
class
atomsmm.integrators.
LimitedSpeedBAOABIntegrator
(stepSize, loops, temperature, frictionConstant, **kwargs)[source]¶ - Parameters
stepSize (unit.Quantity) – The largest time step for numerically integrating the system of equations.
loops (list(int)) – See description in
MultipleTimeScaleIntegrator
.temperature (unit.Quantity) – The temperature to which the configurational sampling should correspond.
frictionConstant (unit.Quantity) – The friction constant \(\gamma\) present in the stochastic equation of motion for per-DOF thermostat variable \(v_2\).
**kwargs (keyword arguments) – The same keyword arguments of class
MultipleTimeScaleIntegrator
apply here.
-
class
atomsmm.integrators.
LimitedSpeedNHLIntegrator
(stepSize, loops, temperature, timeScale, frictionConstant, **kwargs)[source]¶ - Parameters
stepSize (unit.Quantity) – The largest time step for numerically integrating the system of equations.
loops (list(int)) – See description in
MultipleTimeScaleIntegrator
.temperature (unit.Quantity) – The temperature to which the configurational sampling should correspond.
frictionConstant (unit.Quantity) – The friction constant \(\gamma\) present in the stochastic equation of motion for per-DOF thermostat variable \(v_2\).
**kwargs (keyword arguments) – The same keyword arguments of class
MultipleTimeScaleIntegrator
apply here.
-
class
atomsmm.integrators.
LimitedSpeedStochasticIntegrator
(stepSize, loops, temperature, timeScale, frictionConstant, **kwargs)[source]¶ - Parameters
stepSize (unit.Quantity) – The largest time step for numerically integrating the system of equations.
loops (list(int)) – See description in
MultipleTimeScaleIntegrator
.temperature (unit.Quantity) – The temperature to which the configurational sampling should correspond.
frictionConstant (unit.Quantity) – The friction constant \(\gamma\) present in the stochastic equation of motion for per-DOF thermostat variable \(v_2\).
**kwargs (keyword arguments) – The same keyword arguments of class
MultipleTimeScaleIntegrator
apply here.
-
class
atomsmm.integrators.
LimitedSpeedStochasticVelocityIntegrator
(stepSize, loops, temperature, timeScale, frictionConstant, **kwargs)[source]¶ - Parameters
stepSize (unit.Quantity) – The largest time step for numerically integrating the system of equations.
loops (list(int)) – See description in
MultipleTimeScaleIntegrator
.temperature (unit.Quantity) – The temperature to which the configurational sampling should correspond.
frictionConstant (unit.Quantity) – The friction constant \(\gamma\) present in the stochastic equation of motion for per-DOF thermostat variable \(v_2\).
**kwargs (keyword arguments) – The same keyword arguments of class
MultipleTimeScaleIntegrator
apply here.
-
class
atomsmm.integrators.
ExtendedSystemVariable
(name, mass, kT, time_scale, lower_limit=0, upper_limit=1, periodic=False, thermostat='Nose-Hoover', friction_constant=Quantity(value=0.1, unit=/femtosecond))[source]¶ An extended-system variable used for Adiabatic Free Energy Dynamics (AFED).
- Parameters
name (str) – The name of the extended-space variable.
mass (Number or unit.Quantity) – The mass of the extended-space variable.
kT (Number of unit.Quantity) – The temperature of the extended-space variable.
time_scale (Number of unit.Quantity) – The time scale of the thermostat that controls the temperature of this variable.
lower_limit (Number, optional, default=0) – The lower limit of the variable.
upper_limit (Number, optional, default=1) – The upper limit of the variable.
periodic (Bool, optional, default=False) – Whether this variable is subject to periodic boundary conditions. If this is False, then hard, ellastic walls will be considered instead.
-
class
atomsmm.integrators.
AdiabaticDynamicsIntegrator
(custom_integrator, nsteps, variables)[source]¶ This class implements the Adiabatic Free Energy Dynamics (AFED) method.
The equations of motion go as follows:
\[\begin{split}e^{2 n_\mathrm{nsteps} \delta t \mathcal{L}} = & [e^{\frac{1}{2} \delta t \mathbf{F}^t_\lambda \nabla_{\mathbf{p}_\lambda}} e^{\delta t \mathcal{L}_\mathrm{atoms}} e^{\frac{1}{2} \delta t \mathbf{F}^t_\lambda \nabla_{\mathbf{p}_\lambda}} ]^{n_\mathrm{nsteps}} \times \\ & e^{n_\mathrm{nsteps} \delta t \mathbf{p}^t_\lambda \mathbf{M}^{-1}_\lambda \nabla_{\mathbf{\lambda}}} e^{2 n_\mathrm{nsteps} \delta t \mathcal{L}_\mathrm{bath}} e^{n_\mathrm{nsteps} \delta t \mathbf{p}^t_\lambda \mathbf{M}^{-1}_\lambda \nabla_{\mathbf{\lambda}}} \times \\ & [e^{\frac{1}{2} \delta t \mathbf{F}^t_\lambda \nabla_{\mathbf{p}_\lambda}} e^{\delta t \mathcal{L}_\mathrm{atoms}} e^{\frac{1}{2} \delta t \mathbf{F}^t_\lambda \nabla_{\mathbf{p}_\lambda}} ]^{n_\mathrm{nsteps}}\end{split}\]- Parameters
custom_integrator (openmm.CustomIntegrator) – A CustomIntegrator employed to solve the equations of motion of the physical particles, that is, to enact the propagator \(e^{\delta t \mathcal{L}_\mathrm{atoms}}\). The size of an overall AFED time step will be given by \(\Delta t = 2 n_\mathrm{nsteps} \delta t\), where \(\delta t\) is the time step size previously specified for the custom_integrator.
nsteps (int) – The number of consecutive custom_integrator steps executed in the begining of an overall AFED step, and then again in the end.
variables (list(
ExtendedSystemVariable
)) – A list of extended-system variables whose adiabatic dynamics must be taken into account.