DrudeNoseHooverIntegrator

class openmm.openmm.DrudeNoseHooverIntegrator(*args)

This Integrator simulates systems that include Drude particles. It applies two different Nose-Hoover chain thermostats to the different parts of the system. The first is applied to ordinary particles (ones that are not part of a Drude particle pair), as well as to the center of mass of each Drude particle pair. A second thermostat, typically with a much lower temperature, is applied to the relative internal displacement of each pair.

This integrator can optionally set an upper limit on how far any Drude particle is ever allowed to get from its parent particle. This can sometimes help to improve stability. The limit is enforced with a hard wall constraint. By default the limit is set to 0.02 nm.

This Integrator requires the System to include a DrudeForce, which it uses to identify the Drude particles.

__init__(self, temperature, collisionFrequency, drudeTemperature, drudeCollisionFrequency, stepSize, chainLength=3, numMTS=3, numYoshidaSuzuki=7)DrudeNoseHooverIntegrator
__init__(self, other)DrudeNoseHooverIntegrator

Create a DrudeNoseHooverIntegrator.

Parameters
  • temperature (double) – the target temperature for the system (in Kelvin).

  • collisionFrequency (double) – the frequency of the system’s interaction with the heat bath (in inverse picoseconds).

  • drudeTemperature (double) – the target temperature for the Drude particles, relative to their parent atom (in Kelvin).

  • drudeCollisionFrequency (double) – the frequency of the drude particles’ interaction with the heat bath (in inverse picoseconds).

  • stepSize (double) – the step size with which to integrator 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, or 5).

Methods

__init__(-> DrudeNoseHooverIntegrator)

Create a DrudeNoseHooverIntegrator.

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

computeDrudeKineticEnergy(self)

Compute the kinetic energy of the drude particles at the current time.

computeDrudeTemperature(self)

Compute the instantaneous temperature of the Drude system, measured in Kelvin.

computeHeatBathEnergy(self)

Compute the total (potential + kinetic) heat bath energy for all heat baths associated with this integrator, at the current time.

computeSystemTemperature(self)

Compute the instantaneous temperature of the System, measured in Kelvin.

computeTotalKineticEnergy(self)

Compute the kinetic energy of all (real and drude) particles 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.

getIntegrationForceGroups(self)

Get which force groups to use for integration.

getMaxDrudeDistance(self)

Get the maximum distance a Drude particle can ever move from its parent particle, measured in nm.

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.

initialize(self, context)

This will be called by the Context when it is created.

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.

setMaxDrudeDistance(self, distance)

Set the maximum distance a Drude particle can ever move from its parent particle, measured in nm.

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

thisown

The membership flag

property thisown

The membership flag

initialize(self, context)

This will be called by the Context when it is created. It informs the Integrator of what context it will be integrating, and gives it a chance to do any necessary initialization. It will also get called again if the application calls reinitialize() on the Context.

getMaxDrudeDistance(self)double

Get the maximum distance a Drude particle can ever move from its parent particle, measured in nm. This is implemented with a hard wall constraint. If this distance is set to 0 (the default), the hard wall constraint is omitted.

setMaxDrudeDistance(self, distance)

Set the maximum distance a Drude particle can ever move from its parent particle, measured in nm. This is implemented with a hard wall constraint. If this distance is set to 0 (the default), the hard wall constraint is omitted.

computeDrudeKineticEnergy(self)double

Compute the kinetic energy of the drude particles at the current time.

computeTotalKineticEnergy(self)double

Compute the kinetic energy of all (real and drude) particles at the current time.

computeSystemTemperature(self)double

Compute the instantaneous temperature of the System, measured in Kelvin. This is calculated based on the kinetic energy of the ordinary particles (ones not attached to a Drude particle), as well as the center of mass motion of the Drude particle pairs. It does not include the internal motion of the pairs. On average, this should be approximately equal to the value returned by getTemperature().

computeDrudeTemperature(self)double

Compute the instantaneous temperature of the Drude system, measured in Kelvin. This is calculated based on the kinetic energy of the internal motion of Drude pairs and should remain close to the prescribed Drude temperature.

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).

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).

computeHeatBathEnergy(self)double

Compute the total (potential + kinetic) heat bath energy for all heat baths associated with this integrator, at the current time.

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

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.

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.

getNumThermostats(self)int

Get the number of Nose-Hoover chains registered with this integrator.

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

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

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

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

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.

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).

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.

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.

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).

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).

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

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).

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