openmmtools.mcmc.BaseIntegratorMove¶
-
class
openmmtools.mcmc.BaseIntegratorMove(n_steps, context_cache=None, reassign_velocities=False, n_restart_attempts=4)[source]¶ A general MCMC move that applies an integrator.
This class is intended to be inherited by MCMCMoves that need to integrate the system. The child class has to implement the _get_integrator method.
You can decide to override _before_integration() and _after_integration() to execute some code at specific points of the workflow, for example to read data from the Context before the it is destroyed.
Parameters: - n_steps : int
The number of integration steps to take each time the move is applied.
- context_cache : openmmtools.cache.ContextCache, optional
The ContextCache to use for Context creation. If None, the global cache openmmtools.cache.global_context_cache is used (default is None).
- reassign_velocities : bool, optional
If True, the velocities will be reassigned from the Maxwell-Boltzmann distribution at the beginning of the move (default is False).
- restart_attempts : int, optional
When greater than 0, if after the integration there are NaNs in energies, the move will restart. When the integrator has a random component, this may help recovering. On the last attempt, the
Contextis re-initialized in a slower process, but better than the simulation crashing. An IntegratorMoveError is raised after the given number of attempts if there are still NaNs.
Examples
Create a VerletIntegratorMove class.
>>> from openmmtools import testsystems, states >>> from simtk.openmm import VerletIntegrator >>> class VerletMove(BaseIntegratorMove): ... def __init__(self, timestep, n_steps, **kwargs): ... super(VerletMove, self).__init__(n_steps, **kwargs) ... self.timestep = timestep ... def _get_integrator(self, thermodynamic_state): ... return VerletIntegrator(self.timestep) ... def _before_integration(self, context, thermodynamic_state): ... print('Setting velocities') ... context.setVelocitiesToTemperature(thermodynamic_state.temperature) ... def _after_integration(self, context, thermodynamic_state): ... print('Reading statistics') ... >>> alanine = testsystems.AlanineDipeptideVacuum() >>> sampler_state = states.SamplerState(alanine.positions) >>> thermodynamic_state = states.ThermodynamicState(alanine.system, 300*unit.kelvin) >>> move = VerletMove(timestep=1.0*unit.femtosecond, n_steps=2) >>> move.apply(thermodynamic_state, sampler_state) Setting velocities Reading statistics
Attributes: - n_steps : int
- context_cache : openmmtools.cache.ContextCache
- reassign_velocities : bool
- restart_attempts : int or None
Methods
apply(thermodynamic_state, sampler_state)Propagate the state through the integrator. -
__init__(n_steps, context_cache=None, reassign_velocities=False, n_restart_attempts=4)[source]¶ x.__init__(…) initializes x; see help(type(x)) for signature
Methods
__init__(n_steps[, context_cache, …])x.__init__(…) initializes x; see help(type(x)) for signature apply(thermodynamic_state, sampler_state)Propagate the state through the integrator. -
apply(thermodynamic_state, sampler_state)[source]¶ Propagate the state through the integrator.
This updates the SamplerState after the integration. It also logs benchmarking information through the utils.Timer class.
Parameters: - thermodynamic_state : openmmtools.states.ThermodynamicState
The thermodynamic state to use to propagate dynamics.
- sampler_state : openmmtools.states.SamplerState
The sampler state to apply the move to. This is modified.
See also