NoseHooverChain
¶
-
class NoseHooverChain¶
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
G. J. 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.
Public Functions
-
NoseHooverChain(double temperature, double relativeTemperature, double collisionFrequency, double relativeCollisionFrequency, int numDOFs, int chainLength, int numMTS, int numYoshidaSuzuki, int chainID, const std::vector<int> &thermostatedAtoms, const std::vector<std::pair<int, int>> &thermostatedPairs)¶
Create a NoseHooverChain.
- Parameters
temperature – the temperature of the heat bath for absolute motion (in Kelvin)
collisionFrequency – the collision frequency for absolute motion (in 1/ps)
relativeTemperature – 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 – the collision frequency for relative motion(in 1/ps). This is only used if the list of thermostated pairs is not empty.
numDOFs – the number of degrees of freedom in the particles that interact with this chain
chainLength – the length of (number of particles in) this heat bath
numMTS – the number of multi time steps used to propagate this chain
numYoshidaSuzuki – the number of Yoshida Suzuki steps used to propagate this chain (1, 3, 5, or 7).
chainID – 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 – the list of atoms to be handled by this thermostat
thermostatedPairs – 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.
-
inline double getTemperature() const¶
Get the temperature of the heat bath for treating absolute particle motion (in Kelvin).
- Returns
the temperature of the heat bath, measured in Kelvin.
-
inline void setTemperature(double 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 – the temperature of the heat bath (in Kelvin)
-
inline double getRelativeTemperature() const¶
Get the temperature of the heat bath for treating relative particle motion (in Kelvin).
- Returns
the temperature of the heat bath, measured in Kelvin.
-
inline void setRelativeTemperature(double 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 – the temperature of the heat bath for relative motion (in Kelvin)
-
inline double getCollisionFrequency() const¶
Get the collision frequency for treating absolute particle motion (in 1/ps).
- Returns
the collision frequency, measured in 1/ps.
-
inline void setCollisionFrequency(double 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 – the collision frequency (in 1/ps)
-
inline double getRelativeCollisionFrequency() const¶
Get the collision frequency for treating relative particle motion (in 1/ps).
- Returns
the collision frequency, measured in 1/ps.
-
inline void setRelativeCollisionFrequency(double 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 – the collision frequency (in 1/ps)
-
inline int getNumDegreesOfFreedom() const¶
Get the number of degrees of freedom in the particles controled by this heat bath.
- Returns
the number of degrees of freedom.
-
inline void setNumDegreesOfFreedom(int 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 – the number of degrees of freedom.
-
inline int getChainLength() const¶
Get the chain length of this heat bath.
- Returns
the chain length.
-
inline int getNumMultiTimeSteps() const¶
Get the number of steps used in the multi time step propagation.
- Returns
the number of multi time steps.
-
inline int getNumYoshidaSuzukiTimeSteps() const¶
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.
-
inline int getChainID() const¶
Get the chain id used to identify this chain
- Returns
the chain id
-
inline const std::vector<int> &getThermostatedAtoms() const¶
Get the atom ids of all atoms that are thermostated
- Returns
ids of all atoms that are being handled by this thermostat
-
inline void setThermostatedAtoms(const std::vector<int> &atomIDs)¶
Set list of atoms that are handled by this thermostat
- Parameters
atomIDs –
-
inline const std::vector<std::pair<int, int>> &getThermostatedPairs() const¶
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.
-
inline void setThermostatedPairs(const std::vector<std::pair<int, int>> &pairIDs)¶
In case this thermostat handles the kinetic energy of Drude particles set the atom IDs of all parent atoms.
- Parameters
pairIDs – the list of connected pairs to thermostat.
-
std::vector<double> getYoshidaSuzukiWeights() const¶
Get the weights used in the Yoshida Suzuki multi time step decomposition (dimensionless)
- Returns
the weights for the Yoshida-Suzuki integration
-
inline bool usesPeriodicBoundaryConditions() const¶
Returns whether or not this force makes use of periodic boundary conditions.
- Returns
true if force uses PBC and false otherwise
-
NoseHooverChain(double temperature, double relativeTemperature, double collisionFrequency, double relativeCollisionFrequency, int numDOFs, int chainLength, int numMTS, int numYoshidaSuzuki, int chainID, const std::vector<int> &thermostatedAtoms, const std::vector<std::pair<int, int>> &thermostatedPairs)¶