2014-11-25 01:22:05 +08:00
|
|
|
"LAMMPS WWW Site"_lws - "LAMMPS Documentation"_ld - "LAMMPS Commands"_lc :c
|
|
|
|
|
|
|
|
:link(lws,http://lammps.sandia.gov)
|
|
|
|
:link(ld,Manual.html)
|
|
|
|
:link(lc,Section_commands.html#comm)
|
|
|
|
|
|
|
|
:line
|
|
|
|
|
|
|
|
fix gle command :h3
|
|
|
|
|
|
|
|
[Syntax:]
|
|
|
|
|
2015-05-16 00:50:58 +08:00
|
|
|
fix ID id-group gle Ns Tstart Tstop seed Amatrix \[noneq Cmatrix\] \[every stride\] :pre
|
2014-11-25 01:22:05 +08:00
|
|
|
|
|
|
|
ID, group-ID are documented in "fix"_fix.html command :ulb,l
|
|
|
|
gle = style name of this fix command :l
|
|
|
|
Ns = number of additional fictitious momenta :l
|
|
|
|
Tstart, Tstop = temperature ramp during the run :l
|
|
|
|
Amatrix = file to read the drift matrix A from :l
|
|
|
|
seed = random number seed to use for generating noise (positive integer) :l
|
|
|
|
zero or more keyword/value pairs may be appended :l
|
|
|
|
keyword = {noneq} and/or {every}
|
|
|
|
{noneq} Cmatrix = file to read the non-equilibrium covariance matrix from
|
|
|
|
{every} stride = apply the GLE once every time steps. Reduces the accuracy
|
|
|
|
of the integration of the GLE, but has *no effect* on the accuracy of equilibrium
|
|
|
|
sampling. It might change sampling properties when used together with {noneq}.
|
|
|
|
:ule
|
|
|
|
|
|
|
|
[Examples:]
|
|
|
|
|
|
|
|
fix 3 boundary gle 6 300 300 31415 smart.A
|
|
|
|
fix 1 all gle 6 300 300 31415 qt-300k.A noneq qt-300k.C
|
|
|
|
|
|
|
|
[Description:]
|
|
|
|
|
|
|
|
Apply a Generalized Langevin Equation (GLE) thermostat as described
|
|
|
|
in "(Ceriotti)"_#Ceriotti. The formalism allows one to obtain a number
|
|
|
|
of different effects ranging from efficient sampling of all
|
|
|
|
vibrational modes in the system to inexpensive (approximate)
|
|
|
|
modelling of nuclear quantum effects. Contrary to
|
|
|
|
"fix langevin"_fix_langevin.html, this fix performs both
|
|
|
|
thermostatting and evolution of the Hamiltonian equations of motion, so it
|
2014-12-06 06:05:01 +08:00
|
|
|
should not be used together with "fix nve"_fix_nve.html -- at least not
|
|
|
|
on the same atom groups.
|
2014-11-25 01:22:05 +08:00
|
|
|
|
|
|
|
Each degree of freedom in the thermostatted group is supplemented
|
|
|
|
with Ns additional degrees of freedom s, and the equations of motion
|
|
|
|
become
|
|
|
|
|
|
|
|
dq/dt=p/m
|
|
|
|
d(p,s)/dt=(F,0) - A(p,s) + B dW/dt
|
|
|
|
|
|
|
|
where F is the physical force, A is the drift matrix (that generalizes
|
|
|
|
the friction in Langevin dynamics), B is the diffusion term and dW/dt
|
|
|
|
un-correlated Gaussian random forces. The A matrix couples the physical
|
|
|
|
(q,p) dynamics with that of the additional degrees of freedom,
|
|
|
|
and makes it possible to obtain effectively a history-dependent
|
|
|
|
noise and friction kernel.
|
|
|
|
|
|
|
|
The drift matrix should be given as an external file {Afile},
|
|
|
|
as a (Ns+1 x Ns+1) matrix in inverse time units. Matrices that are
|
|
|
|
optimal for a given application and the system of choice can be
|
|
|
|
obtained from "(GLE4MD)"_#GLE4MD.
|
|
|
|
|
|
|
|
Equilibrium sampling a temperature T is obtained by specifiying the
|
|
|
|
target value as the {Tstart} and {Tstop} arguments, so that the diffusion
|
|
|
|
matrix that gives canonical sampling for a given A is computed automatically.
|
|
|
|
However, the GLE framework also allow for non-equilibrium sampling, that
|
|
|
|
can be used for instance to model inexpensively zero-point energy
|
|
|
|
effects "(Ceriotti2)"_#Ceriotti2. This is achieved specifying the
|
|
|
|
{noneq} keyword followed by the name of the file that contains the
|
|
|
|
static covariance matrix for the non-equilibrium dynamics.
|
|
|
|
|
|
|
|
Since integrating GLE dynamics can be costly when used together with
|
|
|
|
simple potentials, one can use the {every} optional keyword to
|
|
|
|
apply the Langevin terms only once every several MD steps, in a
|
|
|
|
multiple time-step fashion. This should be used with care when doing
|
|
|
|
non-equilibrium sampling, but should have no effect on equilibrium
|
|
|
|
averages when using canonical sampling.
|
|
|
|
|
|
|
|
The random number {seed} must be a positive integer. A Marsaglia random
|
|
|
|
number generator is used. Each processor uses the input seed to
|
|
|
|
generate its own unique seed and its own stream of random numbers.
|
|
|
|
Thus the dynamics of the system will not be identical on two runs on
|
|
|
|
different numbers of processors.
|
|
|
|
|
2014-12-06 06:05:01 +08:00
|
|
|
Note also that the Generalized Langevin Dynamics scheme that is
|
|
|
|
implemented by the "fix gld"_fix_gld.html scheme is closely related
|
|
|
|
to the present one. In fact, it should be always possible to cast the
|
|
|
|
Prony series form of the memory kernel used by GLD into an appropriate
|
|
|
|
input matrix for "fix_gle"_fix_gle.html. While the GLE scheme is more
|
|
|
|
general, the form used by "fix gld"_fix_gld.html can be more directly
|
|
|
|
related to the representation of an implicit solvent environment.
|
|
|
|
|
2014-11-25 01:22:05 +08:00
|
|
|
[Restart, fix_modify, output, run start/stop, minimize info:]
|
|
|
|
|
|
|
|
The instantaneous values of the extended variables are written to
|
|
|
|
"binary restart files"_restart.html. Because the state of the random
|
|
|
|
number generator is not saved in restart files, this means you cannot
|
|
|
|
do "exact" restarts with this fix, where the simulation continues on
|
|
|
|
the same as if no restart had taken place. However, in a statistical
|
|
|
|
sense, a restarted simulation should produce the same behavior.
|
|
|
|
Note however that you should use a different seed each time you
|
|
|
|
restart, otherwise the same sequence of random numbers will be used
|
|
|
|
each time, which might lead to stochastic synchronization and
|
|
|
|
subtle artefacts in the sampling.
|
|
|
|
|
|
|
|
This fix can ramp its target temperature over multiple runs, using the
|
|
|
|
{start} and {stop} keywords of the "run"_run.html command. See the
|
|
|
|
"run"_run.html command for details of how to do this.
|
|
|
|
|
|
|
|
The "fix_modify"_fix_modify.html {energy} option is supported by this
|
|
|
|
fix to add the energy change induced by Langevin thermostatting to the
|
|
|
|
system's potential energy as part of "thermodynamic
|
|
|
|
output"_thermo_style.html.
|
|
|
|
|
|
|
|
This fix computes a global scalar which can be accessed by various
|
|
|
|
"output commands"_Section_howto.html#howto_15. The scalar is the
|
|
|
|
cummulative energy change due to this fix. The scalar value
|
|
|
|
calculated by this fix is "extensive".
|
|
|
|
|
|
|
|
[Restrictions:]
|
|
|
|
|
|
|
|
The GLE thermostat in its current implementation should not be used
|
|
|
|
with rigid bodies, SHAKE or RATTLE. It is expected that all the
|
|
|
|
thermostatted degrees of freedom are fully flexible, and the sampled
|
|
|
|
ensemble will not be correct otherwise.
|
|
|
|
|
2014-12-06 06:05:01 +08:00
|
|
|
In order to perform constant-pressure simulations please use
|
|
|
|
"fix press/berendsen"_fix_press_berendsen.html, rather than
|
|
|
|
"fix_npt"_fix_npt.html, to avoid duplicate integration of the
|
|
|
|
equations of motion.
|
|
|
|
|
2014-11-25 01:22:05 +08:00
|
|
|
This fix is part of the USER-MISC package. It is only enabled if LAMMPS
|
|
|
|
was built with that package. See the "Making
|
|
|
|
LAMMPS"_Section_start.html#start_3 section for more info.
|
|
|
|
|
|
|
|
[Related commands:]
|
|
|
|
|
|
|
|
"fix nvt"_fix_nh.html, "fix temp/rescale"_fix_temp_rescale.html, "fix
|
|
|
|
viscous"_fix_viscous.html, "fix nvt"_fix_nh.html, "pair_style
|
|
|
|
dpd/tstat"_pair_dpd.html, "fix_gld"_fix_gld.html
|
|
|
|
|
|
|
|
:line
|
|
|
|
|
|
|
|
:link(Ceriotti)
|
|
|
|
[(Ceriotti)] Ceriotti, Bussi and Parrinello, J Chem Theory Comput 6,
|
|
|
|
1170-80 (2010)
|
|
|
|
|
|
|
|
:link(GLE4MD)
|
|
|
|
[(GLE4MD)] "http://epfl-cosmo.github.io/gle4md/"_http://epfl-cosmo.github.io/gle4md/
|
|
|
|
|
|
|
|
:link(Ceriotti2)
|
|
|
|
[(Ceriotti2)] Ceriotti, Bussi and Parrinello, Phys Rev Lett 103,
|
|
|
|
030603 (2009)
|