NoseHooverChain

class openmm.openmm.NoseHooverChain(*args)

This class defines a chain of Nose-Hoover particles to be used as a heat bath to scale the velocities of a collection of particles subject to thermostating. The heat bath is propagated using the multi time step approach detailed in

    1. Martyna, M. E. Tuckerman, D. J. Tobias and M. L. Klein, Mol. Phys. 87, 1117 (1996).

where the total number of timesteps used to propagate the chain in each step is the number of MTS steps multiplied by the number of terms in the Yoshida-Suzuki decomposition.

Two types of NHC may be created. The first is a simple thermostat that couples with a given subset of the atoms within a system, controling their absolute motion. The second is more elaborate and can thermostat tethered pairs of atoms and in this case two thermostats are created: one that controls the absolute center of mass velocity of each pair and another that controls their motion relative to one another.

__init__(self, temperature, relativeTemperature, collisionFrequency, relativeCollisionFrequency, numDOFs, chainLength, numMTS, numYoshidaSuzuki, chainID, thermostatedAtoms, thermostatedPairs)NoseHooverChain
__init__(self, other)NoseHooverChain

Create a NoseHooverChain.

Parameters
  • temperature (double) – the temperature of the heat bath for absolute motion (in Kelvin)

  • collisionFrequency (double) – the collision frequency for absolute motion (in 1/ps)

  • relativeTemperature (double) – the temperature of the heat bath for relative motion(in Kelvin). This is only used if the list of thermostated pairs is not empty.

  • relativeCollisionFrequency (double) – the collision frequency for relative motion(in 1/ps). This is only used if the list of thermostated pairs is not empty.

  • numDOFs (int) – the number of degrees of freedom in the particles that interact with this chain

  • chainLength (int) – the length of (number of particles in) this heat bath

  • numMTS (int) – the number of multi time steps used to propagate this chain

  • numYoshidaSuzuki (int) – the number of Yoshida Suzuki steps used to propagate this chain (1, 3, 5, or 7).

  • chainID (int) – the chain id used to distinguish this Nose-Hoover chain from others that may be used to control a different set of particles, e.g. for Drude oscillators

  • thermostatedAtoms (vector< int >) – the list of atoms to be handled by this thermostat

  • thermostatedPairs (vector< std::pair< int, int > >) – the list of connected pairs to be thermostated; their absolute center of mass motion will be thermostated independently from their motion relative to one another.

Methods

__init__(-> NoseHooverChain)

Create a NoseHooverChain.

getChainID(self)

Get the chain id used to identify this chain

getChainLength(self)

Get the chain length of this heat bath.

getCollisionFrequency(self)

Get the collision frequency for treating absolute particle motion (in 1/ps).

getNumDegreesOfFreedom(self)

Get the number of degrees of freedom in the particles controled by this heat bath.

getNumMultiTimeSteps(self)

Get the number of steps used in the multi time step propagation.

getNumYoshidaSuzukiTimeSteps(self)

Get the number of steps used in the Yoshida-Suzuki decomposition for multi time step propagation.

getRelativeCollisionFrequency(self)

Get the collision frequency for treating relative particle motion (in 1/ps).

getRelativeTemperature(self)

Get the temperature of the heat bath for treating relative particle motion (in Kelvin).

getTemperature(self)

Get the temperature of the heat bath for treating absolute particle motion (in Kelvin).

getThermostatedAtoms(self)

Get the atom ids of all atoms that are thermostated

getThermostatedPairs(self)

Get the list of any connected pairs to be handled by this thermostat.

getYoshidaSuzukiWeights(self)

Get the weights used in the Yoshida Suzuki multi time step decomposition (dimensionless)

setCollisionFrequency(self, frequency)

Set the collision frequency for treating absolute particle motion.

setNumDegreesOfFreedom(self, numDOF)

Set the number of degrees of freedom in the particles controled by this heat bath.

setRelativeCollisionFrequency(self, frequency)

Set the collision frequency for treating relative particle motion if this thermostat has been set up to handle connected pairs of atoms.

setRelativeTemperature(self, temperature)

Set the temperature of the heat bath for treating relative motion if this thermostat has been set up to treat connected pairs of atoms.

setTemperature(self, temperature)

Set the temperature of the heat bath for treating absolute particle motion.

setThermostatedAtoms(self, atomIDs)

Set list of atoms that are handled by this thermostat

setThermostatedPairs(self, pairIDs)

In case this thermostat handles the kinetic energy of Drude particles set the atom IDs of all parent atoms.

usesPeriodicBoundaryConditions(self)

Returns whether or not this force makes use of periodic boundary conditions.

Attributes

thisown

The membership flag

property thisown

The membership flag

getTemperature(self)double

Get the temperature of the heat bath for treating absolute particle motion (in Kelvin).

Returns

the temperature of the heat bath, measured in Kelvin.

Return type

double

setTemperature(self, temperature)

Set the temperature of the heat bath for treating absolute particle motion. This will affect any new Contexts you create, but not ones that already exist.

Parameters

temperature (double) – the temperature of the heat bath (in Kelvin)

getRelativeTemperature(self)double

Get the temperature of the heat bath for treating relative particle motion (in Kelvin).

Returns

the temperature of the heat bath, measured in Kelvin.

Return type

double

setRelativeTemperature(self, temperature)

Set the temperature of the heat bath for treating relative motion if this thermostat has been set up to treat connected pairs of atoms. This will affect any new Contexts you create, but not ones that already exist.

Parameters

temperature (double) – the temperature of the heat bath for relative motion (in Kelvin)

getCollisionFrequency(self)double

Get the collision frequency for treating absolute particle motion (in 1/ps).

Returns

the collision frequency, measured in 1/ps.

Return type

double

setCollisionFrequency(self, frequency)

Set the collision frequency for treating absolute particle motion. This will affect any new Contexts you create, but not those that already exist.

Parameters

frequency (double) – the collision frequency (in 1/ps)

getRelativeCollisionFrequency(self)double

Get the collision frequency for treating relative particle motion (in 1/ps).

Returns

the collision frequency, measured in 1/ps.

Return type

double

setRelativeCollisionFrequency(self, frequency)

Set the collision frequency for treating relative particle motion if this thermostat has been set up to handle connected pairs of atoms. This will affect any new Contexts you create, but not those that already exist.

Parameters

frequency (double) – the collision frequency (in 1/ps)

getNumDegreesOfFreedom(self)int

Get the number of degrees of freedom in the particles controled by this heat bath.

Returns

the number of degrees of freedom.

Return type

int

setNumDegreesOfFreedom(self, numDOF)

Set the number of degrees of freedom in the particles controled by this heat bath. This will affect any new Contexts you create, but not those that already exist.

Parameters

numDOF (int) – the number of degrees of freedom.

getChainLength(self)int

Get the chain length of this heat bath.

Returns

the chain length.

Return type

int

getNumMultiTimeSteps(self)int

Get the number of steps used in the multi time step propagation.

Returns

the number of multi time steps.

Return type

int

getNumYoshidaSuzukiTimeSteps(self)int

Get the number of steps used in the Yoshida-Suzuki decomposition for multi time step propagation.

Returns

the number of multi time steps in the Yoshida-Suzuki decomposition.

Return type

int

getChainID(self)int

Get the chain id used to identify this chain

Returns

the chain id

Return type

int

getThermostatedAtoms(self)vectori

Get the atom ids of all atoms that are thermostated

Returns

ids of all atoms that are being handled by this thermostat

Return type

vector< int >

setThermostatedAtoms(self, atomIDs)

Set list of atoms that are handled by this thermostat

Parameters

atomIDs (vector< int >) –

getThermostatedPairs(self)vectorpairii

Get the list of any connected pairs to be handled by this thermostat. If this is a regular thermostat, returns an empty vector.

Returns

list of connected pairs.

Return type

vector< std::pair< int, int > >

setThermostatedPairs(self, pairIDs)

In case this thermostat handles the kinetic energy of Drude particles set the atom IDs of all parent atoms.

Parameters

pairIDs (vector< std::pair< int, int > >) – the list of connected pairs to thermostat.

getYoshidaSuzukiWeights(self)vectord

Get the weights used in the Yoshida Suzuki multi time step decomposition (dimensionless)

Returns

the weights for the Yoshida-Suzuki integration

Return type

vector< double >

usesPeriodicBoundaryConditions(self)bool

Returns whether or not this force makes use of periodic boundary conditions.

Returns

true if force uses PBC and false otherwise

Return type

bool