NoseHooverIntegrator¶
- class openmm.openmm.NoseHooverIntegrator(*args)¶
This is an Integrator which simulates a System using one or more Nose Hoover chain thermostats, using the “middle” leapfrog propagation algorithm described in J. Phys. Chem. A 2019, 123, 6056-6079.
- __init__(self, stepSize) → NoseHooverIntegrator¶
- __init__(self, temperature, collisionFrequency, stepSize, chainLength=3, numMTS=3, numYoshidaSuzuki=7) → NoseHooverIntegrator
- __init__(self, other) → NoseHooverIntegrator
Create a NoseHooverIntegrator.
- Parameters
temperature (double) – the target temperature for the system (in Kelvin).
collisionFrequency (double) – the frequency of the interaction with the heat bath (in inverse picoseconds).
stepSize (double) – the step size with which to integrate the system (in picoseconds)
chainLength (int) – the number of beads in the Nose-Hoover chain.
numMTS (int) – the number of step in the multiple time step chain propagation algorithm.
numYoshidaSuzuki (int) – the number of terms in the Yoshida-Suzuki multi time step decomposition used in the chain propagation algorithm (must be 1, 3, 5, or 7).
Methods
__init__
(…)Create a NoseHooverIntegrator.
addSubsystemThermostat
(self, …[, …])Add a Nose-Hoover Chain thermostat to control the temperature of a collection of atoms and/or pairs of connected atoms within the full system.
addThermostat
(self, temperature, …)Add a simple Nose-Hoover Chain thermostat to control the temperature of the full system
computeHeatBathEnergy
(self)Compute the total (potential + kinetic) heat bath energy for all heat baths associated with this integrator, at the current time.
getCollisionFrequency
(self[, chainID])Get the collision frequency for absolute motion of the i-th chain (in 1/picosecond).
getConstraintTolerance
(self)Get the distance tolerance within which constraints are maintained, as a fraction of the constrained distance.
Get which force groups to use for integration.
getMaximumPairDistance
(self)Gets the maximum distance (in nm) that a connected pair may stray from each other.
getNumThermostats
(self)Get the number of Nose-Hoover chains registered with this integrator.
getRelativeCollisionFrequency
(self[, chainID])Get the collision frequency for pairs’ relative motion of the i-th chain (in 1/picosecond).
getRelativeTemperature
(self[, chainID])Get the temperature of the i-th chain for controling pairs’ relative particle motion (in Kelvin).
getStepSize
(self)Get the size of each time step, in picoseconds.
getTemperature
(self[, chainID])Get the temperature of the i-th chain for controling absolute particle motion (in Kelvin).
getThermostat
(self[, chainID])Get the NoseHooverChain thermostat
hasSubsystemThermostats
(self)Return false, if this integrator was set up with the ‘default constructor’ that thermostats the whole system, true otherwise.
setCollisionFrequency
(self, frequency[, chainID])Set the collision frequency for absolute motion of the i-th chain.
setConstraintTolerance
(self, tol)Set the distance tolerance within which constraints are maintained, as a fraction of the constrained distance.
setIntegrationForceGroups
(groups)Set which force groups to use for integration.
setMaximumPairDistance
(self, distance)Sets the maximum distance (in nm) that a connected pair may stray from each other, implemented using a hard wall.
setRelativeCollisionFrequency
(self, frequency)Set the collision frequency for pairs’ relative motion of the i-th chain.
setRelativeTemperature
(self, temperature[, …])set the (relative pair motion) temperature of the i-th chain.
setStepSize
(self, size)Set the size of each time step, in picoseconds.
setTemperature
(self, temperature[, chainID])set the (absolute motion) temperature of the i-th chain.
step
(self, steps)Advance a simulation through time by taking a series of time steps.
Attributes
The membership flag
- property thisown¶
The membership flag
- step(self, steps)¶
Advance a simulation through time by taking a series of time steps.
- Parameters
steps (int) – the number of time steps to take
- addThermostat(self, temperature, collisionFrequency, chainLength, numMTS, numYoshidaSuzuki) → int¶
Add a simple Nose-Hoover Chain thermostat to control the temperature of the full system
- Parameters
temperature (double) – the target temperature for the system.
collisionFrequency (double) – the frequency of the interaction with the heat bath (in 1/ps).
chainLength (int) – the number of beads in the Nose-Hoover chain
numMTS (int) – the number of step in the multiple time step chain propagation algorithm.
numYoshidaSuzuki (int) – the number of terms in the Yoshida-Suzuki multi time step decomposition used in the chain propagation algorithm (must be 1, 3, 5, or 7).
- addSubsystemThermostat(self, thermostatedParticles, thermostatedPairs, temperature, collisionFrequency, relativeTemperature, relativeCollisionFrequency, chainLength=3, numMTS=3, numYoshidaSuzuki=7) → int¶
Add a Nose-Hoover Chain thermostat to control the temperature of a collection of atoms and/or pairs of connected atoms within the full system. A list of atoms defining the atoms to be thermostated is provided and the thermostat will only control members of that list. Additionally a list of pairs of connected atoms may be provided; in this case both the center of mass absolute motion of each pair is controlled as well as their motion relative to each other, which is independently thermostated. If both the list of thermostated particles and thermostated pairs are empty all particles will be thermostated.
- Parameters
thermostatedParticles (vector< int >) – list of particle ids to be thermostated.
thermostatedPairs (vector< std::pair< int, int > >) – a list of pairs of connected atoms whose absolute center of mass motion and motion relative to one another will be independently thermostated.
temperature (double) – the target temperature for each pair’s absolute of center of mass motion.
collisionFrequency (double) – the frequency of the interaction with the heat bath for the pairs’ center of mass motion (in 1/ps).
relativeTemperature (double) – the target temperature for each pair’s relative motion.
relativeCollisionFrequency (double) – the frequency of the interaction with the heat bath for the pairs’ relative motion (in 1/ps).
chainLength (int) – the number of beads in the Nose-Hoover chain.
numMTS (int) – the number of step in the multiple time step chain propagation algorithm.
numYoshidaSuzuki (int) – the number of terms in the Yoshida-Suzuki multi time step decomposition used in the chain propagation algorithm (must be 1, 3, 5, or 7).
- getTemperature(self, chainID=0) → double¶
Get the temperature of the i-th chain for controling absolute particle motion (in Kelvin).
- Parameters
chainID (int) – the index of the Nose-Hoover chain thermostat (default=0).
- Returns
the temperature.
- Return type
double
- setTemperature(self, temperature, chainID=0)¶
set the (absolute motion) temperature of the i-th chain.
- Parameters
temperature (double) – the temperature for the Nose-Hoover chain thermostat (in Kelvin).
chainID (int) – The id of the Nose-Hoover chain thermostat for which the temperature is set (default=0).
- getRelativeTemperature(self, chainID=0) → double¶
Get the temperature of the i-th chain for controling pairs’ relative particle motion (in Kelvin).
- Parameters
chainID (int) – the index of the Nose-Hoover chain thermostat (default=0).
- Returns
the temperature.
- Return type
double
- setRelativeTemperature(self, temperature, chainID=0)¶
set the (relative pair motion) temperature of the i-th chain.
- Parameters
temperature (double) – the temperature for the Nose-Hoover chain thermostat (in Kelvin).
chainID (int) – The id of the Nose-Hoover chain thermostat for which the temperature is set (default=0).
- getCollisionFrequency(self, chainID=0) → double¶
Get the collision frequency for absolute motion of the i-th chain (in 1/picosecond).
- Parameters
chainID (int) – the index of the Nose-Hoover chain thermostat (default=0).
- Returns
the collision frequency.
- Return type
double
- setCollisionFrequency(self, frequency, chainID=0)¶
Set the collision frequency for absolute motion of the i-th chain.
- Parameters
frequency (double) – the collision frequency in picosecond.
chainID (int) – the index of the Nose-Hoover chain thermostat (default=0).
- getRelativeCollisionFrequency(self, chainID=0) → double¶
Get the collision frequency for pairs’ relative motion of the i-th chain (in 1/picosecond).
- Parameters
chainID (int) – the index of the Nose-Hoover chain thermostat (default=0).
- Returns
the collision frequency.
- Return type
double
- setRelativeCollisionFrequency(self, frequency, chainID=0)¶
Set the collision frequency for pairs’ relative motion of the i-th chain.
- Parameters
frequency (double) – the collision frequency in picosecond.
chainID (int) – the index of the Nose-Hoover chain thermostat (default=0).
- computeHeatBathEnergy(self) → double¶
Compute the total (potential + kinetic) heat bath energy for all heat baths associated with this integrator, at the current time.
- getNumThermostats(self) → int¶
Get the number of Nose-Hoover chains registered with this integrator.
- getThermostat(self, chainID=0) → NoseHooverChain¶
Get the NoseHooverChain thermostat
- Parameters
chainID (int) – the index of the Nose-Hoover chain thermostat (default=0).
- hasSubsystemThermostats(self) → bool¶
Return false, if this integrator was set up with the ‘default constructor’ that thermostats the whole system, true otherwise. Required for serialization.
- getMaximumPairDistance(self) → double¶
Gets the maximum distance (in nm) that a connected pair may stray from each other. If zero, there are no constraints on the intra-pair separation.
- setMaximumPairDistance(self, distance)¶
Sets the maximum distance (in nm) that a connected pair may stray from each other, implemented using a hard wall. If set to zero, the hard wall constraint is omited and the pairs are free to be separated by any distance.
- getConstraintTolerance(self) → double¶
Get the distance tolerance within which constraints are maintained, as a fraction of the constrained distance.
- getIntegrationForceGroups(self) → int¶
Get which force groups to use for integration. By default, all force groups are included. This is interpreted as a set of bit flags: the forces from group i will be included if (groups&(1<<i)) != 0.
- getStepSize(self) → double¶
Get the size of each time step, in picoseconds. If this integrator uses variable time steps, the size of the most recent step is returned.
- Returns
the step size, measured in ps
- Return type
double
- setConstraintTolerance(self, tol)¶
Set the distance tolerance within which constraints are maintained, as a fraction of the constrained distance.
- setIntegrationForceGroups(groups)¶
Set which force groups to use for integration. By default, all force groups are included.
- Parameters
groups (set or int) – a set of indices for which force groups to include when integrating the equations of motion. Alternatively, the groups can be passed as a single unsigned integer interpreted as a bitmask, in which case group i will be included if (groups&(1<<i)) != 0.
- setStepSize(self, size)¶
Set the size of each time step, in picoseconds. If this integrator uses variable time steps, the effect of calling this method is undefined, and it may simply be ignored.
- Parameters
size (double) – the step size, measured in ps