NoseHooverIntegrator

class OpenMM::NoseHooverIntegrator

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.

Methods

NoseHooverIntegrator

Create a NoseHooverIntegrator.

NoseHooverIntegrator

Create a NoseHooverIntegrator.

~NoseHooverIntegrator

step

Advance a simulation through time by taking a series of time steps.

addThermostat

Add a simple Nose-Hoover Chain thermostat to control the temperature of the full system

addSubsystemThermostat

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.

getTemperature

Get the temperature of the i-th chain for controling absolute particle motion (in Kelvin).

setTemperature

set the (absolute motion) temperature of the i-th chain.

getRelativeTemperature

Get the temperature of the i-th chain for controling pairs’ relative particle motion (in Kelvin).

setRelativeTemperature

set the (relative pair motion) temperature of the i-th chain.

getCollisionFrequency

Get the collision frequency for absolute motion of the i-th chain (in 1/picosecond).

setCollisionFrequency

Set the collision frequency for absolute motion of the i-th chain.

getRelativeCollisionFrequency

Get the collision frequency for pairs’ relative motion of the i-th chain (in 1/picosecond).

setRelativeCollisionFrequency

Set the collision frequency for pairs’ relative motion of the i-th chain.

computeHeatBathEnergy

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

getNumThermostats

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

getThermostat

Get the NoseHooverChain thermostat

hasSubsystemThermostats

Return false, if this integrator was set up with the ‘default constructor’ that thermostats the whole system, true otherwise.

getMaximumPairDistance

Gets the maximum distance (in nm) that a connected pair may stray from each other.

setMaximumPairDistance

Sets the maximum distance (in nm) that a connected pair may stray from each other, implemented using a hard wall.

getAllThermostatedIndividualParticles

Get a list of all individual atoms (i.e.

getAllThermostatedPairs

Get a list of all connected Drude-like pairs, and their target relative temperature, in the system.

NoseHooverIntegrator(double stepSize)

Create a NoseHooverIntegrator(). This version creates a bare leapfrog integrator with no thermostats; any thermostats should be added by calling addThermostat.

Parameters

  • stepSize – the step size with which to integrate the system (in picoseconds)

NoseHooverIntegrator(double temperature, double collisionFrequency, double stepSize, int chainLength = 3, int numMTS = 3, int numYoshidaSuzuki = 7)

Create a NoseHooverIntegrator().

Parameters

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

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

  • stepSize – the step size with which to integrate the system (in picoseconds)

  • chainLength – the number of beads in the Nose-Hoover chain.

  • numMTS – the number of step in the multiple time step chain propagation algorithm.

  • numYoshidaSuzuki – 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).

~NoseHooverIntegrator()
void step(int steps)

Advance a simulation through time by taking a series of time steps.

Parameters

  • steps – the number of time steps to take

int addThermostat(double temperature, double collisionFrequency, int chainLength, int numMTS, int numYoshidaSuzuki)

Add a simple Nose-Hoover Chain thermostat to control the temperature of the full system

Parameters

  • temperature – the target temperature for the system.

  • collisionFrequency – the frequency of the interaction with the heat bath (in 1/ps).

  • chainLength – the number of beads in the Nose-Hoover chain

  • numMTS – the number of step in the multiple time step chain propagation algorithm.

  • numYoshidaSuzuki – 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).

int addSubsystemThermostat(const std::vector<int> &thermostatedParticles, const std::vector<std::pair<int, int>> &thermostatedPairs, double temperature, double collisionFrequency, double relativeTemperature, double relativeCollisionFrequency, int chainLength = 3, int numMTS = 3, int numYoshidaSuzuki = 7)

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 – list of particle ids to be thermostated.

  • thermostatedPairs – a list of pairs of connected atoms whose absolute center of mass motion and motion relative to one another will be independently thermostated.

  • temperature – the target temperature for each pair’s absolute of center of mass motion.

  • collisionFrequency – the frequency of the interaction with the heat bath for the pairs’ center of mass motion (in 1/ps).

  • relativeTemperature – the target temperature for each pair’s relative motion.

  • relativeCollisionFrequency – the frequency of the interaction with the heat bath for the pairs’ relative motion (in 1/ps).

  • chainLength – the number of beads in the Nose-Hoover chain.

  • numMTS – the number of step in the multiple time step chain propagation algorithm.

  • numYoshidaSuzuki – 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).

double getTemperature(int chainID = 0) const

Get the temperature of the i-th chain for controling absolute particle motion (in Kelvin).

Parameters

  • chainID – the index of the Nose-Hoover chain thermostat (default=0).

Returns

the temperature.

void setTemperature(double temperature, int chainID = 0)

set the (absolute motion) temperature of the i-th chain.

Parameters

  • temperature – the temperature for the Nose-Hoover chain thermostat (in Kelvin).

  • chainID – The id of the Nose-Hoover chain thermostat for which the temperature is set (default=0).

double getRelativeTemperature(int chainID = 0) const

Get the temperature of the i-th chain for controling pairs’ relative particle motion (in Kelvin).

Parameters

  • chainID – the index of the Nose-Hoover chain thermostat (default=0).

Returns

the temperature.

void setRelativeTemperature(double temperature, int chainID = 0)

set the (relative pair motion) temperature of the i-th chain.

Parameters

  • temperature – the temperature for the Nose-Hoover chain thermostat (in Kelvin).

  • chainID – The id of the Nose-Hoover chain thermostat for which the temperature is set (default=0).

double getCollisionFrequency(int chainID = 0) const

Get the collision frequency for absolute motion of the i-th chain (in 1/picosecond).

Parameters

  • chainID – the index of the Nose-Hoover chain thermostat (default=0).

Returns

the collision frequency.

void setCollisionFrequency(double frequency, int chainID = 0)

Set the collision frequency for absolute motion of the i-th chain.

Parameters

  • frequency – the collision frequency in picosecond.

  • chainID – the index of the Nose-Hoover chain thermostat (default=0).

double getRelativeCollisionFrequency(int chainID = 0) const

Get the collision frequency for pairs’ relative motion of the i-th chain (in 1/picosecond).

Parameters

  • chainID – the index of the Nose-Hoover chain thermostat (default=0).

Returns

the collision frequency.

void setRelativeCollisionFrequency(double frequency, int chainID = 0)

Set the collision frequency for pairs’ relative motion of the i-th chain.

Parameters

  • frequency – the collision frequency in picosecond.

  • chainID – the index of the Nose-Hoover chain thermostat (default=0).

double computeHeatBathEnergy()

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

int getNumThermostats() const

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

const NoseHooverChain &getThermostat(int chainID = 0) const

Get the NoseHooverChain thermostat

Parameters

  • chainID – the index of the Nose-Hoover chain thermostat (default=0).

bool hasSubsystemThermostats() const

Return false, if this integrator was set up with the ‘default constructor’ that thermostats the whole system, true otherwise. Required for serialization.

double getMaximumPairDistance() const

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.

void setMaximumPairDistance(double 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.

const std::vector<int> &getAllThermostatedIndividualParticles() const

Get a list of all individual atoms (i.e. not involved in a connected Drude-like pair) in the system.

const std::vector<std::tuple<int, int, double>> &getAllThermostatedPairs() const

Get a list of all connected Drude-like pairs, and their target relative temperature, in the system.