git-svn-id: svn://svn.icms.temple.edu/lammps-ro/trunk@16053 f3b2605a-c512-4ea7-a41b...
[lammps.git] / doc / src / fix_nh.txt
blob0938509b45dd02a16fa9a3338fb7096722fa3ebb
1 <"LAMMPS WWW Site"_lws - "LAMMPS Documentation"_ld - "LAMMPS Commands"_lc :c
3 :link(lws,http://lammps.sandia.gov)
4 :link(ld,Manual.html)
5 :link(lc,Section_commands.html#comm)
7 :line
9 fix nvt command :h3
10 fix nvt/intel command :h3
11 fix nvt/kk command :h3
12 fix nvt/omp command :h3
13 fix npt command :h3
14 fix npt/intel command :h3
15 fix npt/kk command :h3
16 fix npt/omp command :h3
17 fix nph command :h3
18 fix nph/kk command :h3
19 fix nph/omp command :h3
21 [Syntax:]
23 fix ID group-ID style_name keyword value ... :pre
25 ID, group-ID are documented in "fix"_fix.html command :ulb,l
26 style_name = {nvt} or {npt} or {nph} :l
27 one or more keyword/value pairs may be appended :l
28 keyword = {temp} or {iso} or {aniso} or {tri} or {x} or {y} or {z} or {xy} or {yz} or {xz} or {couple} or {tchain} or {pchain} or {mtk} or {tloop} or {ploop} or {nreset} or {drag} or {dilate} or {scalexy} or {scaleyz} or {scalexz} or {flip} or {fixedpoint} or {update}
29   {temp} values = Tstart Tstop Tdamp
30     Tstart,Tstop = external temperature at start/end of run
31     Tdamp = temperature damping parameter (time units)
32   {iso} or {aniso} or {tri} values = Pstart Pstop Pdamp
33     Pstart,Pstop = scalar external pressure at start/end of run (pressure units)
34     Pdamp = pressure damping parameter (time units)
35   {x} or {y} or {z} or {xy} or {yz} or {xz} values = Pstart Pstop Pdamp
36     Pstart,Pstop = external stress tensor component at start/end of run (pressure units)
37     Pdamp = stress damping parameter (time units)
38   {couple} = {none} or {xyz} or {xy} or {yz} or {xz}
39   {tchain} value = N
40     N = length of thermostat chain (1 = single thermostat)
41   {pchain} values = N
42     N length of thermostat chain on barostat (0 = no thermostat)
43   {mtk} value = {yes} or {no} = add in MTK adjustment term or not
44   {tloop} value = M
45     M = number of sub-cycles to perform on thermostat
46   {ploop} value = M
47     M = number of sub-cycles to perform on barostat thermostat
48   {nreset} value = reset reference cell every this many timesteps
49   {drag} value = Df
50     Df = drag factor added to barostat/thermostat (0.0 = no drag)
51   {dilate} value = dilate-group-ID
52     dilate-group-ID = only dilate atoms in this group due to barostat volume changes
53   {scalexy} value = {yes} or {no} = scale xy with ly
54   {scaleyz} value = {yes} or {no} = scale yz with lz
55   {scalexz} value = {yes} or {no} = scale xz with lz
56   {flip} value = {yes} or {no} = allow or disallow box flips when it becomes highly skewed
57   {fixedpoint} values = x y z
58     x,y,z = perform barostat dilation/contraction around this point (distance units)
59   {update} value = {dipole} or {dipole/dlm}
60     dipole = update dipole orientation (only for sphere variants)
61     dipole/dlm = use DLM integrator to update dipole orientation (only for sphere variants) :pre
62 :ule
64 [Examples:]
66 fix 1 all nvt temp 300.0 300.0 100.0
67 fix 1 water npt temp 300.0 300.0 100.0 iso 0.0 0.0 1000.0
68 fix 2 jello npt temp 300.0 300.0 100.0 tri 5.0 5.0 1000.0
69 fix 2 ice nph x 1.0 1.0 0.5 y 2.0 2.0 0.5 z 3.0 3.0 0.5 yz 0.1 0.1 0.5 xz 0.2 0.2 0.5 xy 0.3 0.3 0.5 nreset 1000 :pre
71 [Description:]
73 These commands perform time integration on Nose-Hoover style
74 non-Hamiltonian equations of motion which are designed to generate
75 positions and velocities sampled from the canonical (nvt),
76 isothermal-isobaric (npt), and isenthalpic (nph) ensembles.  This
77 updates the position and velocity for atoms in the group each
78 timestep.
80 The thermostatting and barostatting is achieved by adding some dynamic
81 variables which are coupled to the particle velocities
82 (thermostatting) and simulation domain dimensions (barostatting).  In
83 addition to basic thermostatting and barostatting, these fixes can
84 also create a chain of thermostats coupled to the particle thermostat,
85 and another chain of thermostats coupled to the barostat
86 variables. The barostat can be coupled to the overall box volume, or
87 to individual dimensions, including the {xy}, {xz} and {yz} tilt
88 dimensions. The external pressure of the barostat can be specified as
89 either a scalar pressure (isobaric ensemble) or as components of a
90 symmetric stress tensor (constant stress ensemble).  When used
91 correctly, the time-averaged temperature and stress tensor of the
92 particles will match the target values specified by Tstart/Tstop and
93 Pstart/Pstop.
95 The equations of motion used are those of Shinoda et al in
96 "(Shinoda)"_#nh-Shinoda, which combine the hydrostatic equations of
97 Martyna, Tobias and Klein in "(Martyna)"_#nh-Martyna with the strain
98 energy proposed by Parrinello and Rahman in
99 "(Parrinello)"_#nh-Parrinello.  The time integration schemes closely
100 follow the time-reversible measure-preserving Verlet and rRESPA
101 integrators derived by Tuckerman et al in "(Tuckerman)"_#nh-Tuckerman.
103 :line
105 The thermostat parameters for fix styles {nvt} and {npt} is specified
106 using the {temp} keyword.  Other thermostat-related keywords are
107 {tchain}, {tloop} and {drag}, which are discussed below.
109 The thermostat is applied to only the translational degrees of freedom
110 for the particles.  The translational degrees of freedom can also have
111 a bias velocity removed before thermostatting takes place; see the
112 description below.  The desired temperature at each timestep is a
113 ramped value during the run from {Tstart} to {Tstop}.  The {Tdamp}
114 parameter is specified in time units and determines how rapidly the
115 temperature is relaxed.  For example, a value of 10.0 means to relax
116 the temperature in a timespan of (roughly) 10 time units (e.g. tau or
117 fmsec or psec - see the "units"_units.html command).  The atoms in the
118 fix group are the only ones whose velocities and positions are updated
119 by the velocity/position update portion of the integration.
121 NOTE: A Nose-Hoover thermostat will not work well for arbitrary values
122 of {Tdamp}.  If {Tdamp} is too small, the temperature can fluctuate
123 wildly; if it is too large, the temperature will take a very long time
124 to equilibrate.  A good choice for many models is a {Tdamp} of around
125 100 timesteps.  Note that this is NOT the same as 100 time units for
126 most "units"_units.html settings.
128 :line
130 The barostat parameters for fix styles {npt} and {nph} is specified
131 using one or more of the {iso}, {aniso}, {tri}, {x}, {y}, {z}, {xy},
132 {xz}, {yz}, and {couple} keywords.  These keywords give you the
133 ability to specify all 6 components of an external stress tensor, and
134 to couple various of these components together so that the dimensions
135 they represent are varied together during a constant-pressure
136 simulation.
138 Other barostat-related keywords are {pchain}, {mtk}, {ploop},
139 {nreset}, {drag}, and {dilate}, which are discussed below.
141 Orthogonal simulation boxes have 3 adjustable dimensions (x,y,z).
142 Triclinic (non-orthogonal) simulation boxes have 6 adjustable
143 dimensions (x,y,z,xy,xz,yz).  The "create_box"_create_box.html, "read
144 data"_read_data.html, and "read_restart"_read_restart.html commands
145 specify whether the simulation box is orthogonal or non-orthogonal
146 (triclinic) and explain the meaning of the xy,xz,yz tilt factors.
148 The target pressures for each of the 6 components of the stress tensor
149 can be specified independently via the {x}, {y}, {z}, {xy}, {xz}, {yz}
150 keywords, which correspond to the 6 simulation box dimensions.  For
151 each component, the external pressure or tensor component at each
152 timestep is a ramped value during the run from {Pstart} to {Pstop}.
153 If a target pressure is specified for a component, then the
154 corresponding box dimension will change during a simulation.  For
155 example, if the {y} keyword is used, the y-box length will change.  If
156 the {xy} keyword is used, the xy tilt factor will change.  A box
157 dimension will not change if that component is not specified, although
158 you have the option to change that dimension via the "fix
159 deform"_fix_deform.html command.
161 Note that in order to use the {xy}, {xz}, or {yz} keywords, the
162 simulation box must be triclinic, even if its initial tilt factors are
163 0.0.
165 For all barostat keywords, the {Pdamp} parameter operates like the
166 {Tdamp} parameter, determining the time scale on which pressure is
167 relaxed.  For example, a value of 10.0 means to relax the pressure in
168 a timespan of (roughly) 10 time units (e.g. tau or fmsec or psec - see
169 the "units"_units.html command).
171 NOTE: A Nose-Hoover barostat will not work well for arbitrary values
172 of {Pdamp}.  If {Pdamp} is too small, the pressure and volume can
173 fluctuate wildly; if it is too large, the pressure will take a very
174 long time to equilibrate.  A good choice for many models is a {Pdamp}
175 of around 1000 timesteps.  However, note that {Pdamp} is specified in
176 time units, and that timesteps are NOT the same as time units for most
177 "units"_units.html settings.
179 Regardless of what atoms are in the fix group (the only atoms which
180 are time integrated), a global pressure or stress tensor is computed
181 for all atoms.  Similarly, when the size of the simulation box is
182 changed, all atoms are re-scaled to new positions, unless the keyword
183 {dilate} is specified with a {dilate-group-ID} for a group that
184 represents a subset of the atoms.  This can be useful, for example, to
185 leave the coordinates of atoms in a solid substrate unchanged and
186 controlling the pressure of a surrounding fluid.  This option should
187 be used with care, since it can be unphysical to dilate some atoms and
188 not others, because it can introduce large, instantaneous
189 displacements between a pair of atoms (one dilated, one not) that are
190 far from the dilation origin.  Also note that for atoms not in the fix
191 group, a separate time integration fix like "fix nve"_fix_nve.html or
192 "fix nvt"_fix_nh.html can be used on them, independent of whether they
193 are dilated or not.
195 :line
197 The {couple} keyword allows two or three of the diagonal components of
198 the pressure tensor to be "coupled" together.  The value specified
199 with the keyword determines which are coupled.  For example, {xz}
200 means the {Pxx} and {Pzz} components of the stress tensor are coupled.
201 {Xyz} means all 3 diagonal components are coupled.  Coupling means two
202 things: the instantaneous stress will be computed as an average of the
203 corresponding diagonal components, and the coupled box dimensions will
204 be changed together in lockstep, meaning coupled dimensions will be
205 dilated or contracted by the same percentage every timestep.  The
206 {Pstart}, {Pstop}, {Pdamp} parameters for any coupled dimensions must
207 be identical.  {Couple xyz} can be used for a 2d simulation; the {z}
208 dimension is simply ignored.
210 :line
212 The {iso}, {aniso}, and {tri} keywords are simply shortcuts that are
213 equivalent to specifying several other keywords together.
215 The keyword {iso} means couple all 3 diagonal components together when
216 pressure is computed (hydrostatic pressure), and dilate/contract the
217 dimensions together.  Using "iso Pstart Pstop Pdamp" is the same as
218 specifying these 4 keywords:
220 x Pstart Pstop Pdamp
221 y Pstart Pstop Pdamp
222 z Pstart Pstop Pdamp
223 couple xyz :pre
225 The keyword {aniso} means {x}, {y}, and {z} dimensions are controlled
226 independently using the {Pxx}, {Pyy}, and {Pzz} components of the
227 stress tensor as the driving forces, and the specified scalar external
228 pressure.  Using "aniso Pstart Pstop Pdamp" is the same as specifying
229 these 4 keywords:
231 x Pstart Pstop Pdamp
232 y Pstart Pstop Pdamp
233 z Pstart Pstop Pdamp
234 couple none :pre
236 The keyword {tri} means {x}, {y}, {z}, {xy}, {xz}, and {yz} dimensions
237 are controlled independently using their individual stress components
238 as the driving forces, and the specified scalar pressure as the
239 external normal stress.  Using "tri Pstart Pstop Pdamp" is the same as
240 specifying these 7 keywords:
242 x Pstart Pstop Pdamp
243 y Pstart Pstop Pdamp
244 z Pstart Pstop Pdamp
245 xy 0.0 0.0 Pdamp
246 yz 0.0 0.0 Pdamp
247 xz 0.0 0.0 Pdamp
248 couple none :pre
250 :line
252 In some cases (e.g. for solids) the pressure (volume) and/or
253 temperature of the system can oscillate undesirably when a Nose/Hoover
254 barostat and thermostat is applied.  The optional {drag} keyword will
255 damp these oscillations, although it alters the Nose/Hoover equations.
256 A value of 0.0 (no drag) leaves the Nose/Hoover formalism unchanged.
257 A non-zero value adds a drag term; the larger the value specified, the
258 greater the damping effect.  Performing a short run and monitoring the
259 pressure and temperature is the best way to determine if the drag term
260 is working.  Typically a value between 0.2 to 2.0 is sufficient to
261 damp oscillations after a few periods. Note that use of the drag
262 keyword will interfere with energy conservation and will also change
263 the distribution of positions and velocities so that they do not
264 correspond to the nominal NVT, NPT, or NPH ensembles.
266 An alternative way to control initial oscillations is to use chain
267 thermostats. The keyword {tchain} determines the number of thermostats
268 in the particle thermostat. A value of 1 corresponds to the original
269 Nose-Hoover thermostat. The keyword {pchain} specifies the number of
270 thermostats in the chain thermostatting the barostat degrees of
271 freedom. A value of 0 corresponds to no thermostatting of the
272 barostat variables.
274 The {mtk} keyword controls whether or not the correction terms due to
275 Martyna, Tuckerman, and Klein are included in the equations of motion
276 "(Martyna)"_#nh-Martyna.  Specifying {no} reproduces the original
277 Hoover barostat, whose volume probability distribution function
278 differs from the true NPT and NPH ensembles by a factor of 1/V.  Hence
279 using {yes} is more correct, but in many cases the difference is
280 negligible.
282 The keyword {tloop} can be used to improve the accuracy of integration
283 scheme at little extra cost.  The initial and final updates of the
284 thermostat variables are broken up into {tloop} substeps, each of
285 length {dt}/{tloop}. This corresponds to using a first-order
286 Suzuki-Yoshida scheme "(Tuckerman)"_#nh-Tuckerman.  The keyword {ploop}
287 does the same thing for the barostat thermostat.
289 The keyword {nreset} controls how often the reference dimensions used
290 to define the strain energy are reset.  If this keyword is not used,
291 or is given a value of zero, then the reference dimensions are set to
292 those of the initial simulation domain and are never changed. If the
293 simulation domain changes significantly during the simulation, then
294 the final average pressure tensor will differ significantly from the
295 specified values of the external stress tensor.  A value of {nstep}
296 means that every {nstep} timesteps, the reference dimensions are set
297 to those of the current simulation domain.
299 The {scaleyz}, {scalexz}, and {scalexy} keywords control whether or
300 not the corresponding tilt factors are scaled with the associated box
301 dimensions when barostatting triclinic periodic cells.  The default
302 values {yes} will turn on scaling, which corresponds to adjusting the
303 linear dimensions of the cell while preserving its shape.  Choosing
304 {no} ensures that the tilt factors are not scaled with the box
305 dimensions. See below for restrictions and default values in different
306 situations. In older versions of LAMMPS, scaling of tilt factors was
307 not performed. The old behavior can be recovered by setting all three
308 scale keywords to {no}.
310 The {flip} keyword allows the tilt factors for a triclinic box to
311 exceed half the distance of the parallel box length, as discussed
312 below.  If the {flip} value is set to {yes}, the bound is enforced by
313 flipping the box when it is exceeded.  If the {flip} value is set to
314 {no}, the tilt will continue to change without flipping.  Note that if
315 applied stress induces large deformations (e.g. in a liquid), this
316 means the box shape can tilt dramatically and LAMMPS will run less
317 efficiently, due to the large volume of communication needed to
318 acquire ghost atoms around a processor's irregular-shaped sub-domain.
319 For extreme values of tilt, LAMMPS may also lose atoms and generate an
320 error.
322 The {fixedpoint} keyword specifies the fixed point for barostat volume
323 changes. By default, it is the center of the box.  Whatever point is
324 chosen will not move during the simulation.  For example, if the lower
325 periodic boundaries pass through (0,0,0), and this point is provided
326 to {fixedpoint}, then the lower periodic boundaries will remain at
327 (0,0,0), while the upper periodic boundaries will move twice as
328 far. In all cases, the particle trajectories are unaffected by the
329 chosen value, except for a time-dependent constant translation of
330 positions.
332 If the {update} keyword is used with the {dipole} value, then the
333 orientation of the dipole moment of each particle is also updated
334 during the time integration.  This option should be used for models
335 where a dipole moment is assigned to finite-size particles,
336 e.g. spheroids via use of the "atom_style hybrid sphere
337 dipole"_atom_style.html command.
339 The default dipole orientation integrator can be changed to the
340 Dullweber-Leimkuhler-McLachlan integration scheme
341 "(Dullweber)"_#nh-Dullweber when using {update} with the value
342 {dipole/dlm}. This integrator is symplectic and time-reversible,
343 giving better energy conservation and allows slightly longer timesteps
344 at only a small additional computational cost.
346 :line
348 NOTE: Using a barostat coupled to tilt dimensions {xy}, {xz}, {yz} can
349 sometimes result in arbitrarily large values of the tilt dimensions,
350 i.e. a dramatically deformed simulation box.  LAMMPS allows the tilt
351 factors to grow a small amount beyond the normal limit of half the box
352 length (0.6 times the box length), and then performs a box "flip" to
353 an equivalent periodic cell.  See the discussion of the {flip} keyword
354 above, to allow this bound to be exceeded, if desired.
356 The flip operation is described in more detail in the doc page for
357 "fix deform"_fix_deform.html.  Both the barostat dynamics and the atom
358 trajectories are unaffected by this operation.  However, if a tilt
359 factor is incremented by a large amount (1.5 times the box length) on
360 a single timestep, LAMMPS can not accomodate this event and will
361 terminate the simulation with an error. This error typically indicates
362 that there is something badly wrong with how the simulation was
363 constructed, such as specifying values of {Pstart} that are too far
364 from the current stress value, or specifying a timestep that is too
365 large. Triclinic barostatting should be used with care. This also is
366 true for other barostat styles, although they tend to be more
367 forgiving of insults. In particular, it is important to recognize that
368 equilibrium liquids can not support a shear stress and that
369 equilibrium solids can not support shear stresses that exceed the
370 yield stress.
372 One exception to this rule is if the 1st dimension in the tilt factor
373 (x for xy) is non-periodic.  In that case, the limits on the tilt
374 factor are not enforced, since flipping the box in that dimension does
375 not change the atom positions due to non-periodicity.  In this mode,
376 if you tilt the system to extreme angles, the simulation will simply
377 become inefficient due to the highly skewed simulation box.
379 NOTE: Unlike the "fix temp/berendsen"_fix_temp_berendsen.html command
380 which performs thermostatting but NO time integration, these fixes
381 perform thermostatting/barostatting AND time integration.  Thus you
382 should not use any other time integration fix, such as "fix
383 nve"_fix_nve.html on atoms to which this fix is applied.  Likewise,
384 fix nvt and fix npt should not normally be used on atoms that also
385 have their temperature controlled by another fix - e.g. by "fix
386 langevin"_fix_nh.html or "fix temp/rescale"_fix_temp_rescale.html
387 commands.
389 See "this howto section"_Section_howto.html#howto_16 of the manual for
390 a discussion of different ways to compute temperature and perform
391 thermostatting and barostatting.
393 :line
395 These fixes compute a temperature and pressure each timestep.  To do
396 this, the fix creates its own computes of style "temp" and "pressure",
397 as if one of these two sets of commands had been issued:
399 compute fix-ID_temp group-ID temp
400 compute fix-ID_press group-ID pressure fix-ID_temp :pre
402 compute fix-ID_temp all temp
403 compute fix-ID_press all pressure fix-ID_temp :pre
405 See the "compute temp"_compute_temp.html and "compute
406 pressure"_compute_pressure.html commands for details.  Note that the
407 IDs of the new computes are the fix-ID + underscore + "temp" or fix_ID
408 + underscore + "press".  For fix nvt, the group for the new computes
409 is the same as the fix group.  For fix nph and fix npt, the group for
410 the new computes is "all" since pressure is computed for the entire
411 system.
413 Note that these are NOT the computes used by thermodynamic output (see
414 the "thermo_style"_thermo_style.html command) with ID = {thermo_temp}
415 and {thermo_press}.  This means you can change the attributes of this
416 fix's temperature or pressure via the
417 "compute_modify"_compute_modify.html command or print this temperature
418 or pressure during thermodynamic output via the "thermo_style
419 custom"_thermo_style.html command using the appropriate compute-ID.
420 It also means that changing attributes of {thermo_temp} or
421 {thermo_press} will have no effect on this fix.
423 Like other fixes that perform thermostatting, fix nvt and fix npt can
424 be used with "compute commands"_compute.html that calculate a
425 temperature after removing a "bias" from the atom velocities.
426 E.g. removing the center-of-mass velocity from a group of atoms or
427 only calculating temperature on the x-component of velocity or only
428 calculating temperature for atoms in a geometric region.  This is not
429 done by default, but only if the "fix_modify"_fix_modify.html command
430 is used to assign a temperature compute to this fix that includes such
431 a bias term.  See the doc pages for individual "compute
432 commands"_compute.html to determine which ones include a bias.  In
433 this case, the thermostat works in the following manner: the current
434 temperature is calculated taking the bias into account, bias is
435 removed from each atom, thermostatting is performed on the remaining
436 thermal degrees of freedom, and the bias is added back in.
438 :line
440 These fixes can be used with either the {verlet} or {respa}
441 "integrators"_run_style.html. When using one of the barostat fixes
442 with {respa}, LAMMPS uses an integrator constructed
443 according to the following factorization of the Liouville propagator
444 (for two rRESPA levels):
446 :c,image(Eqs/fix_nh1.jpg)
448 This factorization differs somewhat from that of Tuckerman et al, in
449 that the barostat is only updated at the outermost rRESPA level,
450 whereas Tuckerman's factorization requires splitting the pressure into
451 pieces corresponding to the forces computed at each rRESPA level. In
452 theory, the latter method will exhibit better numerical stability. In
453 practice, because Pdamp is normally chosen to be a large multiple of
454 the outermost rRESPA timestep, the barostat dynamics are not the
455 limiting factor for numerical stability. Both factorizations are
456 time-reversible and can be shown to preserve the phase space measure
457 of the underlying non-Hamiltonian equations of motion.
459 NOTE: This implementation has been shown to conserve linear momentum
460 up to machine precision under NVT dynamics. Under NPT dynamics,
461 for a system with zero initial total linear momentum, the total
462 momentum fluctuates close to zero. It may occasionally undergo brief
463 excursions to non-negligible values, before returning close to zero.
464 Over long simulations, this has the effect of causing the center-of-mass
465 to undergo a slow random walk. This can be mitigated by resetting
466 the momentum at infrequent intervals using the
467 "fix momentum"_fix_momentum.html command.
469 NOTE: This implementation has been shown to conserve linear momentum
470 up to machine precision under NVT dynamics. Under NPT dynamics,
471 for a system with zero initial total linear momentum, the total
472 momentum fluctuates close to zero. It may occasionally undergo brief
473 excursions to non-negligible values, before returning close to zero.
474 Over long simulations, this has the effect of causing the center-of-mass
475 to undergo a slow random walk. This can be mitigated by resetting
476 the momentum at infrequent intervals using the
477 "fix momentum"_fix_momentum.html command.
479 :line
481 The fix npt and fix nph commands can be used with rigid bodies or
482 mixtures of rigid bodies and non-rigid particles (e.g. solvent).  But
483 there are also "fix rigid/npt"_fix_rigid.html and "fix
484 rigid/nph"_fix_rigid.html commands, which are typically a more natural
485 choice.  See the doc page for those commands for more discussion of
486 the various ways to do this.
488 :line
490 Styles with a {gpu}, {intel}, {kk}, {omp}, or {opt} suffix are
491 functionally the same as the corresponding style without the suffix.
492 They have been optimized to run faster, depending on your available
493 hardware, as discussed in "Section 5"_Section_accelerate.html
494 of the manual.  The accelerated styles take the same arguments and
495 should produce the same results, except for round-off and precision
496 issues.
498 These accelerated styles are part of the GPU, USER-INTEL, KOKKOS,
499 USER-OMP and OPT packages, respectively.  They are only enabled if
500 LAMMPS was built with those packages.  See the "Making
501 LAMMPS"_Section_start.html#start_3 section for more info.
503 You can specify the accelerated styles explicitly in your input script
504 by including their suffix, or you can use the "-suffix command-line
505 switch"_Section_start.html#start_7 when you invoke LAMMPS, or you can
506 use the "suffix"_suffix.html command in your input script.
508 See "Section 5"_Section_accelerate.html of the manual for
509 more instructions on how to use the accelerated styles effectively.
511 :line
513 [Restart, fix_modify, output, run start/stop, minimize info:]
515 These fixes writes the state of all the thermostat and barostat
516 variables to "binary restart files"_restart.html.  See the
517 "read_restart"_read_restart.html command for info on how to re-specify
518 a fix in an input script that reads a restart file, so that the
519 operation of the fix continues in an uninterrupted fashion.
521 The "fix_modify"_fix_modify.html {temp} and {press} options are
522 supported by these fixes.  You can use them to assign a
523 "compute"_compute.html you have defined to this fix which will be used
524 in its thermostatting or barostatting procedure, as described above.
525 If you do this, note that the kinetic energy derived from the compute
526 temperature should be consistent with the virial term computed using
527 all atoms for the pressure.  LAMMPS will warn you if you choose to
528 compute temperature on a subset of atoms.
530 NOTE: If both the {temp} and {press} keywords are used in a single
531 thermo_modify command (or in two separate commands), then the order in
532 which the keywords are specified is important.  Note that a "pressure
533 compute"_compute_pressure.html defines its own temperature compute as
534 an argument when it is specified.  The {temp} keyword will override
535 this (for the pressure compute being used by fix npt), but only if the
536 {temp} keyword comes after the {press} keyword.  If the {temp} keyword
537 comes before the {press} keyword, then the new pressure compute
538 specified by the {press} keyword will be unaffected by the {temp}
539 setting.
541 The "fix_modify"_fix_modify.html {energy} option is supported by these
542 fixes to add the energy change induced by Nose/Hoover thermostatting
543 and barostatting to the system's potential energy as part of
544 "thermodynamic output"_thermo_style.html.
546 These fixes compute a global scalar and a global vector of quantities,
547 which can be accessed by various "output
548 commands"_Section_howto.html#howto_15.  The scalar value calculated by
549 these fixes is "extensive"; the vector values are "intensive".
551 The scalar is the cumulative energy change due to the fix.
553 The vector stores internal Nose/Hoover thermostat and barostat
554 variables.  The number and meaning of the vector values depends on
555 which fix is used and the settings for keywords {tchain} and {pchain},
556 which specify the number of Nose/Hoover chains for the thermostat and
557 barostat.  If no thermostatting is done, then {tchain} is 0.  If no
558 barostatting is done, then {pchain} is 0.  In the following list,
559 "ndof" is 0, 1, 3, or 6, and is the number of degrees of freedom in
560 the barostat.  Its value is 0 if no barostat is used, else its value
561 is 6 if any off-diagonal stress tensor component is barostatted, else
562 its value is 1 if {couple xyz} is used or {couple xy} for a 2d
563 simulation, otherwise its value is 3.
565 The order of values in the global vector and their meaning is as
566 follows.  The notation means there are tchain values for eta, followed
567 by tchain for eta_dot, followed by ndof for omega, etc:
569 eta\[tchain\] = particle thermostat displacements (unitless)
570 eta_dot\[tchain\] = particle thermostat velocities (1/time units)
571 omega\[ndof\] = barostat displacements (unitless)
572 omega_dot\[ndof\] = barostat velocities (1/time units)
573 etap\[pchain\] = barostat thermostat displacements (unitless)
574 etap_dot\[pchain\] = barostat thermostat velocities (1/time units)
575 PE_eta\[tchain\] = potential energy of each particle thermostat displacement (energy units)
576 KE_eta_dot\[tchain\] = kinetic energy of each particle thermostat velocity (energy units)
577 PE_omega\[ndof\] = potential energy of each barostat displacement (energy units)
578 KE_omega_dot\[ndof\] = kinetic energy of each barostat velocity (energy units)
579 PE_etap\[pchain\] = potential energy of each barostat thermostat displacement (energy units)
580 KE_etap_dot\[pchain\] = kinetic energy of each barostat thermostat velocity (energy units)
581 PE_strain\[1\] = scalar strain energy (energy units) :ul
583 These fixes can ramp their external temperature and pressure over
584 multiple runs, using the {start} and {stop} keywords of the
585 "run"_run.html command.  See the "run"_run.html command for details of
586 how to do this.
588 These fixes are not invoked during "energy
589 minimization"_minimize.html.
591 :line
593 [Restrictions:]
595 {X}, {y}, {z} cannot be barostatted if the associated dimension is not
596 periodic.  {Xy}, {xz}, and {yz} can only be barostatted if the
597 simulation domain is triclinic and the 2nd dimension in the keyword
598 ({y} dimension in {xy}) is periodic.  {Z}, {xz}, and {yz}, cannot be
599 barostatted for 2D simulations.  The "create_box"_create_box.html,
600 "read data"_read_data.html, and "read_restart"_read_restart.html
601 commands specify whether the simulation box is orthogonal or
602 non-orthogonal (triclinic) and explain the meaning of the xy,xz,yz
603 tilt factors.
605 For the {temp} keyword, the final Tstop cannot be 0.0 since it would
606 make the external T = 0.0 at some timestep during the simulation which
607 is not allowed in the Nose/Hoover formulation.
609 The {scaleyz yes} and {scalexz yes} keyword/value pairs can not be used
610 for 2D simulations. {scaleyz yes}, {scalexz yes}, and {scalexy yes} options
611 can only be used if the 2nd dimension in the keyword is periodic,
612 and if the tilt factor is not coupled to the barostat via keywords
613 {tri}, {yz}, {xz}, and {xy}.
615 These fixes can be used with dynamic groups as defined by the
616 "group"_group.html command.  Likewise they can be used with groups to
617 which atoms are added or deleted over time, e.g. a deposition
618 simulation.  However, the conservation properties of the thermostat
619 and barostat are defined for systems with a static set of atoms.  You
620 may observe odd behavior if the atoms in a group vary dramatically
621 over time or the atom count becomes very small.
623 [Related commands:]
625 "fix nve"_fix_nve.html, "fix_modify"_fix_modify.html,
626 "run_style"_run_style.html
628 [Default:]
630 The keyword defaults are tchain = 3, pchain = 3, mtk = yes, tloop =
631 ploop = 1, nreset = 0, drag = 0.0, dilate = all, couple = none,
632 scaleyz = scalexz = scalexy = yes if periodic in 2nd dimension and
633 not coupled to barostat, otherwise no.
635 :line
637 :link(nh-Martyna)
638 [(Martyna)] Martyna, Tobias and Klein, J Chem Phys, 101, 4177 (1994).
640 :link(nh-Parrinello)
641 [(Parrinello)] Parrinello and Rahman, J Appl Phys, 52, 7182 (1981).
643 :link(nh-Tuckerman)
644 [(Tuckerman)] Tuckerman, Alejandre, Lopez-Rendon, Jochim, and
645 Martyna, J Phys A: Math Gen, 39, 5629 (2006).
647 :link(nh-Shinoda)
648 [(Shinoda)] Shinoda, Shiga, and Mikami, Phys Rev B, 69, 134103 (2004).
650 :link(nh-Dullweber)
651 [(Dullweber)] Dullweber, Leimkuhler and McLachlan, J Chem Phys, 107,
652 5840 (1997).