md.methods¶
Overview
Brownian dynamics. |
|
Langevin dynamics. |
|
NPT Integration via MTK barostat-thermostat. |
|
NVE Integration via Velocity-Verlet |
|
NVT Integration via the Nosé-Hoover thermostat. |
Details
-
class
hoomd.md.methods.
Brownian
(filter, kT, seed, alpha=None)¶ Brownian dynamics.
- Parameters
filter (
hoomd.filter.ParticleFilter
) – Subset of particles to apply this method to.kT (
hoomd.variant.Variant
orfloat
) – Temperature of the simulation (in energy units).seed (
int
) – Random seed to use for generating \(\vec{F}_\mathrm{R}\).alpha (
float
) – When set, use \(\alpha d_i\) for the drag coefficient where \(d_i\) is particle diameter. Defaults to None.
Brownian
integrates particles forward in time according to the overdamped Langevin equations of motion, sometimes called Brownian dynamics, or the diffusive limit.\[ \begin{align}\begin{aligned}\frac{d\vec{x}}{dt} = \frac{\vec{F}_\mathrm{C} + \vec{F}_\mathrm{R}}{\gamma}\\\langle \vec{F}_\mathrm{R} \rangle = 0\\\langle |\vec{F}_\mathrm{R}|^2 \rangle = 2 d k T \gamma / \delta t\\\langle \vec{v}(t) \rangle = 0\\\langle |\vec{v}(t)|^2 \rangle = d k T / m\end{aligned}\end{align} \]where \(\vec{F}_\mathrm{C}\) is the force on the particle from all potentials and constraint forces, \(\gamma\) is the drag coefficient, \(\vec{F}_\mathrm{R}\) is a uniform random force, \(\vec{v}\) is the particle’s velocity, and \(d\) is the dimensionality of the system. The magnitude of the random force is chosen via the fluctuation-dissipation theorem to be consistent with the specified drag and temperature, \(T\). When \(kT=0\), the random force \(\vec{F}_\mathrm{R}=0\).
Brownian
generates random numbers by hashing together the particle tag, user seed, and current time step index. See C. L. Phillips et. al. 2011 for more information.Attention
Change the seed if you reset the simulation time step to 0. If you keep the same seed, the simulation will continue with the same sequence of random numbers used previously and may cause unphysical correlations.
For MPI runs: all ranks other than 0 ignore the seed input and use the value of rank 0.
Brownian
uses the integrator from I. Snook, The Langevin and Generalised Langevin Approach to the Dynamics of Atomic, Polymeric and Colloidal Systems , 2007, section 6.2.5, with the exception that \(\vec{F}_\mathrm{R}\) is drawn from a uniform random number distribution.In Brownian dynamics, particle velocities are completely decoupled from positions. At each time step,
Brownian
draws a new velocity distribution consistent with the current set temperature so thathoomd.compute.thermo
will report appropriate temperatures and pressures if logged or needed by other commands.Brownian dynamics neglects the acceleration term in the Langevin equation. This assumption is valid when overdamped: \(\frac{m}{\gamma} \ll \delta t\). Use
Langevin
if your system is not overdamped.You can specify \(\gamma\) in two ways:
Specify \(\alpha\) which scales the particle diameter to \(\gamma = \alpha d_i\). The units of \(\alpha\) are mass / distance / time.
After the method object is created, specify the attribute
gamma
andgamma_r
for rotational damping or random torque to assign them directly, with independent values for each particle type in the system.
Examples:
brownian = hoomd.md.methods.Brownian(filter=hoomd.filter.All(), kT=0.2, seed=1, alpha=1.0) integrator = hoomd.md.Integrator(dt=0.001, methods=[brownian], forces=[lj])
Examples of using
gamma
prgamma_r
on drag coefficient:brownian = hoomd.md.methods.Brownian(filter=hoomd.filter.All(), kT=0.2, seed=1) brownian.gamma.default = 2.0 brownian.gamma_r.default = [1.0, 2.0, 3.0]
-
filter
¶ Subset of particles to apply this method to.
-
kT
¶ Temperature of the simulation (in energy units).
-
alpha
¶ When set, use \(\alpha d_i\) for the drag coefficient where \(d_i\) is particle diameter. Defaults to None.
- Type
-
class
hoomd.md.methods.
Langevin
(filter, kT, seed, alpha=None, tally_reservoir_energy=False)¶ Langevin dynamics.
- Parameters
filter (
hoomd.filter.ParticleFilter
) – Subset of particles to apply this method to.kT (
hoomd.variant.Variant
orfloat
) – Temperature of the simulation (in energy units).seed (
int
) – Random seed to use for generating \(\vec{F}_\mathrm{R}\).alpha (
float
) – When set, use \(\alpha d_i\) for the drag coefficient where \(d_i\) is particle diameter. Defaults to None.tally_reservoir_energy (
bool
) – If true, the energy exchange between the thermal reservoir and the particles is tracked. Total energy conservation can then be monitored by addinglangevin_reservoir_energy_groupname
to the logged quantities. Defaults to False.
Translational degrees of freedom
Langevin
integrates particles forward in time according to the Langevin equations of motion:\[ \begin{align}\begin{aligned}m \frac{d\vec{v}}{dt} = \vec{F}_\mathrm{C} - \gamma \cdot \vec{v} + \vec{F}_\mathrm{R}\\\langle \vec{F}_\mathrm{R} \rangle = 0\\\langle |\vec{F}_\mathrm{R}|^2 \rangle = 2 d kT \gamma / \delta t\end{aligned}\end{align} \]where \(\vec{F}_\mathrm{C}\) is the force on the particle from all potentials and constraint forces, \(\gamma\) is the drag coefficient, \(\vec{v}\) is the particle’s velocity, \(\vec{F}_\mathrm{R}\) is a uniform random force, and \(d\) is the dimensionality of the system (2 or 3). The magnitude of the random force is chosen via the fluctuation-dissipation theorem to be consistent with the specified drag and temperature, \(T\). When \(kT=0\), the random force \(\vec{F}_\mathrm{R}=0\).
Langevin
generates random numbers by hashing together the particle tag, user seed, and current time step index. See C. L. Phillips et. al. 2011 for more information.Attention
Change the seed if you reset the simulation time step to 0. If you keep the same seed, the simulation will continue with the same sequence of random numbers used previously and may cause unphysical correlations.
For MPI runs: all ranks other than 0 ignore the seed input and use the value of rank 0.
Langevin dynamics includes the acceleration term in the Langevin equation and is useful for gently thermalizing systems using a small gamma. This assumption is valid when underdamped: \(\frac{m}{\gamma} \gg \delta t\). Use
Brownian
if your system is not underdamped.Langevin
uses the same integrator asNVE
with the additional force term \(- \gamma \cdot \vec{v} + \vec{F}_\mathrm{R}\). The random force \(\vec{F}_\mathrm{R}\) is drawn from a uniform random number distribution.You can specify \(\gamma\) in two ways:
Specify \(\alpha\) which scales the particle diameter to \(\gamma = \alpha d_i\). The units of \(\alpha\) are mass / distance / time.
After the method object is created, specify the attribute
gamma
andgamma_r
for rotational damping or random torque to assign them directly, with independent values for each particle type in the system.
Warning
When restarting a simulation, the energy of the reservoir will be reset to zero.
Examples:
langevin = hoomd.md.methods.Langevin(filter=hoomd.filter.All(), kT=0.2, seed=1, alpha=1.0) integrator = hoomd.md.Integrator(dt=0.001, methods=[langevin], forces=[lj])
Examples of using
gamma
orgamma_r
on drag coefficient:langevin = hoomd.md.methods.Langevin(filter=hoomd.filter.All(), kT=0.2, seed=1) langevin.gamma.default = 2.0 langevin.gamma_r.default = [1.0,2.0,3.0]
-
filter
¶ Subset of particles to apply this method to.
-
kT
¶ Temperature of the simulation (in energy units).
-
alpha
¶ When set, use \(\alpha d_i\) for the drag coefficient where \(d_i\) is particle diameter. Defaults to None.
- Type
-
class
hoomd.md.methods.
NPT
(filter, kT, tau, S, tauS, couple, box_dof=[True, True, True, False, False, False], rescale_all=False, gamma=0.0)¶ NPT Integration via MTK barostat-thermostat.
- Parameters
filter (
hoomd.filter.ParticleFilter
) – Subset of particles on which to apply this method.kT (
hoomd.variant.Variant
orfloat
) – Temperature set point for the thermostat. (in energy units).tau (
float
) – Coupling constant for the thermostat (in time units).S (
list
[hoomd.variant.Variant
] orfloat
) – Stress components set point for the barostat (in pressure units). In Voigt notation: \([S_{xx}, S_{yy}, S_{zz}, S_{yz}, S_{xz}, S_{xy}]\). In case of isotropic pressure P (\([p, p, p, 0, 0, 0]\)), useS = p
.tauS (
float
) – Coupling constant for the barostat (in time units).couple (
str
) – Couplings of diagonal elements of the stress tensor, can be “none”, “xy”, “xz”,”yz”, or “all”, default to “all”.box_dof (
list
[bool
]) – Box degrees of freedom with six boolean elements corresponding to x, y, z, xy, xz, yz, each. Default to [True,True,True,False,False,False]). If turned on to True, rescale corresponding lengths or tilt factors and components of particle coordinates and velocities.rescale_all (
bool
) – if True, rescale all particles, not only those in the group, Default to False.gamma (
float
) – Dimensionless damping factor for the box degrees of freedom, Default to 0.
NPT
performs constant pressure, constant temperature simulations, allowing for a fully deformable simulation box.The integration method is based on the rigorous Martyna-Tobias-Klein equations of motion for NPT. For optimal stability, the update equations leave the phase-space measure invariant and are manifestly time-reversible.
By default,
NPT
performs integration in a cubic box under hydrostatic pressure by simultaneously rescaling the lengths Lx, Ly and Lz of the simulation box.NPT
can also perform more advanced integration modes. The integration mode is specified by a set of couplings and by specifying the box degrees of freedom that are put under barostat control.Couplings define which diagonal elements of the pressure tensor \(P_{\alpha,\beta}\) should be averaged over, so that the corresponding box lengths are rescaled by the same amount.
Valid couplings are:
none (all box lengths are updated independently)
xy (Lx and Ly are coupled)
xz (Lx and Lz are coupled)
yz (Ly and Lz are coupled)
all (Lx and Ly (and Lz if 3D) are coupled)
The default coupling is all, i.e. the ratios between all box lengths stay constant.
Degrees of freedom of the box specify which lengths and tilt factors of the box should be updated, and how particle coordinates and velocities should be rescaled. The
box_dof
tuple controls the way the box is rescaled and updated. The first three elementsbox_dof[:3]
controls whether the x, y, and z box lengths are rescaled and updated, respectively. The last three entriesbox_dof[3:]
control the rescaling or the tilt factors xy, xz, and yz. All options also appropriately rescale particle coordinates and velocities.By default, the x, y, and z degrees of freedom are updated.
[True,True,True,False,False,False]
Note
If any of the diagonal x, y, z degrees of freedom is not being integrated, pressure tensor components along that direction are not considered for the remaining degrees of freedom.
For example:
Specifying all couplings and x, y, and z degrees of freedom amounts to cubic symmetry (default)
Specifying xy couplings and x, y, and z degrees of freedom amounts to tetragonal symmetry.
Specifying no couplings and all degrees of freedom amounts to a fully deformable triclinic unit cell
For the MTK equations of motion, see:
Glaser et. al (2013), unpublished
\(\tau\) is related to the Nosé mass \(Q\) by
\[\tau = \sqrt{\frac{Q}{g k T_0}}\]where \(g\) is the number of degrees of freedom, and \(k T_0\) is the set point (kT above).
The
NPT
equations of motion include a translational thermostat (with momentum \(\xi\) and position \(\eta\)), a rotational thermostat (with momentum \(\xi_{\mathrm{rot}}\) and position \(\eta_\mathrm{rot}\)), and a barostat tensor \(\nu_{\mathrm{ij}}\). Access these quantities usingtranslational_thermostat_dof
,rotational_thermostat_dof
, andbarostat_dof
.Note
Coupling constant for barostat
tauS
should be set within appropriate range for pressure and volume to fluctuate in reasonable rate and equilibrate. Too smalltauS
can cause abrupt fluctuation, whereas too largetauS
would take long time to equilibrate. In most of systems, recommended value fortauS
is1000 * dt
, wheredt
is the length of the time step.Examples:
npt = hoomd.md.methods.NPT(filter=hoomd.filter.All(), tau=1.0, kT=0.65, tauS = 1.2, S=2.0) # orthorhombic symmetry npt = hoomd.md.methods.NPT(filter=hoomd.filter.All(), tau=1.0, kT=0.65, tauS = 1.2, S=2.0, couple="none") # tetragonal symmetry npt = hoomd.md.methods.NPT(filter=hoomd.filter.All(), tau=1.0, kT=0.65, tauS = 1.2, S=2.0, couple="xy") # triclinic symmetry npt = hoomd.md.methods.NPT(filter=hoomd.filter.All(), tau=1.0, kT=0.65, tauS = 1.2, S=2.0, couple="none", rescale_all=True) integrator = hoomd.md.Integrator(dt=0.005, methods=[npt], forces=[lj])
-
filter
¶ Subset of particles on which to apply this method.
-
kT
¶ Temperature set point for the thermostat. (in energy units).
-
S
¶ Stress components set point for the barostat (in pressure units). In Voigt notation, \([S_{xx}, S_{yy}, S_{zz}, S_{yz}, S_{xz}, S_{xy}]\). Stress can be reset after method object is created. For example, An isoropic pressure can be set by
npt.S = 4.
- Type
List[hoomd.variant.Variant]
-
couple
¶ Couplings of diagonal elements of the stress tensor, can be “none”, “xy”, “xz”,”yz”, or “all”.
- Type
-
box_dof
¶ Box degrees of freedom with six boolean elements corresponding to x, y, z, xy, xz, yz, each.
- Type
List[bool]
-
translational_thermostat_dof
¶ Additional degrees of freedom for the translational thermostat (\(\xi\), \(\eta\))
-
rotational_thermostat_dof
¶ Additional degrees of freedom for the rotational thermostat (\(\xi_\mathrm{rot}\), \(\eta_\mathrm{rot}\))
-
barostat_dof
¶ Additional degrees of freedom for the barostat (\(\nu_{xx}\), \(\nu_{xy}\), \(\nu_{xz}\), \(\nu_{yy}\), \(\nu_{yz}\), \(\nu_{zz}\))
-
thermalize_thermostat_and_barostat_dof
(seed)¶ Set the thermostat and barostat momenta to random values.
- Parameters
seed (int) – Random number seed
thermalize_thermostat_and_barostat_dof
sets a random value for the momentum \(\xi\) and the barostat \(\nu_{\mathrm{ij}}\). WhenIntegrator.aniso
isTrue
, it also sets a random value for the rotational thermostat momentum \(\xi_{\mathrm{rot}}\). Callthermalize_thermostat_and_barostat_dof
to set a new random state for the thermostat and barostat.Important
You must call
Simulation.run
beforethermalize_thermostat_and_barostat_dof
. Callrun(steps=0)
to prepare a newly createdSimulation
.See also
Note
The seed for the pseudorandom number stream includes the simulation timestep and the provided seed.
-
class
hoomd.md.methods.
NVE
(filter, limit=None)¶ NVE Integration via Velocity-Verlet
- Parameters
filter (
hoomd.filter.ParticleFilter
) – Subset of particles on which to apply this method.limit (None or
float
) – Enforce that no particle moves more than a distance of a limit in a single time step. Defaults to None
NVE
performs constant volume, constant energy simulations using the standard Velocity-Verlet method. For poor initial conditions that include overlapping atoms, a limit can be specified to the movement a particle is allowed to make in one time step. After a few thousand time steps with the limit set, the system should be in a safe state to continue with unconstrained integration.Examples:
nve = hoomd.md.methods.NVE(filter=hoomd.filter.All()) integrator = hoomd.md.Integrator(dt=0.005, methods=[nve], forces=[lj])
-
filter
¶ Subset of particles on which to apply this method.
-
class
hoomd.md.methods.
NVT
(filter, kT, tau)¶ NVT Integration via the Nosé-Hoover thermostat.
- Parameters
filter (
hoomd.filter.ParticleFilter
) – Subset of particles on which to apply this method.kT (
hoomd.variant.Variant
orfloat
) – Temperature set point for the Nosé-Hoover thermostat. (in energy units).tau (
float
) – Coupling constant for the Nosé-Hoover thermostat. (in time units).
NVT
performs constant volume, constant temperature simulations using the Nosé-Hoover thermostat, using the MTK equations described in Refs. G. J. Martyna, D. J. Tobias, M. L. Klein 1994 and J. Cao, G. J. Martyna 1996.\(\tau\) is related to the Nosé mass \(Q\) by
\[\tau = \sqrt{\frac{Q}{g k T_0}}\]where \(g\) is the number of degrees of freedom, and \(k T_0\) is the set point (kT above).
The
NVT
equations of motion include a translational thermostat (with momentum \(\xi\) and position \(\eta\)) and a rotational thermostat (with momentum \(\xi_{\mathrm{rot}}\) and position \(\eta_\mathrm{rot}\)). Access these quantities usingtranslational_thermostat_dof
androtational_thermostat_dof
.Note
Coupling constant
tau
in Nosé-Hoover thermostat should be set within reasonable range to avoid abrupt fluctuation in temperature in case of smalltau
, also to avoid long time to equilibrate in case of largetau
. Recommended value for most of systems is100 * dt
, wheredt
is the length of the time step.Examples:
nvt=hoomd.md.methods.NVT(filter=hoomd.filter.All(), kT=1.0, tau=0.5) integrator = hoomd.md.Integrator(dt=0.005, methods=[nvt], forces=[lj])
-
filter
¶ Subset of particles on which to apply this method.
-
kT
¶ Temperature set point for the Nosé-Hoover thermostat. (in energy units).
-
translational_thermostat_dof
¶ Additional degrees of freedom for the translational thermostat (\(\xi\), \(\eta\))
-
rotational_thermostat_dof
¶ Additional degrees of freedom for the rotational thermostat (\(\xi_\mathrm{rot}\), \(\eta_\mathrm{rot}\))
-
thermalize_thermostat_dof
(seed)¶ Set the thermostat momenta to random values.
- Parameters
seed (int) – Random number seed
thermalize_extra_dof
sets a random value for the momentum \(\xi\). WhenIntegrator.aniso
isTrue
, it also sets a random value for the rotational thermostat momentum \(\xi_{\mathrm{rot}}\). Callthermalize_extra_dof
to set a new random state for the thermostat.Important
You must call
Simulation.run
beforethermalize_extra_dof
. Callrun(steps=0)
to prepare a newly createdSimulation
.See also
Note
The seed for the pseudorandom number stream includes the simulation timestep and the provided seed.