Make Cygwin build work
[gromacs.git] / docs / user-guide / mdp-options.rst
blob9281ac4e0e0c7b0730de6705f293479d5d2fa02d
1 .. README
2    See the "run control" section for a working example of the
3    syntax to use when making .mdp entries, with and without detailed
4    documentation for values those entries might take. Everything can
5    be cross-referenced, see the examples there.
7 .. todo:: Make more cross-references.
9 Molecular dynamics parameters (.mdp options)
10 ============================================
12 .. _mdp-general:
14 General information
15 -------------------
17 Default values are given in parentheses, or listed first among
18 choices. The first option in the list is always the default
19 option. Units are given in square brackets. The difference between a
20 dash and an underscore is ignored.
22 A :ref:`sample mdp file <mdp>` is available. This should be
23 appropriate to start a normal simulation. Edit it to suit your
24 specific needs and desires.
27 Preprocessing
28 ^^^^^^^^^^^^^
30 .. mdp:: include
32    directories to include in your topology. Format:
33    ``-I/home/john/mylib -I../otherlib``
35 .. mdp:: define
37    defines to pass to the preprocessor, default is no defines. You can
38    use any defines to control options in your customized topology
39    files. Options that act on existing :ref:`top` file mechanisms
40    include
42       ``-DFLEXIBLE`` will use flexible water instead of rigid water
43       into your topology, this can be useful for normal mode analysis.
45       ``-DPOSRES`` will trigger the inclusion of ``posre.itp`` into
46       your topology, used for implementing position restraints.
49 Run control
50 ^^^^^^^^^^^
52 .. mdp:: integrator
54    (Despite the name, this list includes algorithms that are not
55    actually integrators over time. :mdp-value:`integrator=steep` and
56    all entries following it are in this category)
58    .. mdp-value:: md
60       A leap-frog algorithm for integrating Newton's equations of motion.
62    .. mdp-value:: md-vv
64       A velocity Verlet algorithm for integrating Newton's equations
65       of motion.  For constant NVE simulations started from
66       corresponding points in the same trajectory, the trajectories
67       are analytically, but not binary, identical to the
68       :mdp-value:`integrator=md` leap-frog integrator. The kinetic
69       energy, which is determined from the whole step velocities and
70       is therefore slightly too high. The advantage of this integrator
71       is more accurate, reversible Nose-Hoover and Parrinello-Rahman
72       coupling integration based on Trotter expansion, as well as
73       (slightly too small) full step velocity output. This all comes
74       at the cost off extra computation, especially with constraints
75       and extra communication in parallel. Note that for nearly all
76       production simulations the :mdp-value:`integrator=md` integrator
77       is accurate enough.
79    .. mdp-value:: md-vv-avek
81       A velocity Verlet algorithm identical to
82       :mdp-value:`integrator=md-vv`, except that the kinetic energy is
83       determined as the average of the two half step kinetic energies
84       as in the :mdp-value:`integrator=md` integrator, and this thus
85       more accurate.  With Nose-Hoover and/or Parrinello-Rahman
86       coupling this comes with a slight increase in computational
87       cost.
89    .. mdp-value:: sd
91       An accurate and efficient leap-frog stochastic dynamics
92       integrator. With constraints, coordinates needs to be
93       constrained twice per integration step. Depending on the
94       computational cost of the force calculation, this can take a
95       significant part of the simulation time. The temperature for one
96       or more groups of atoms (:mdp:`tc-grps`) is set with
97       :mdp:`ref-t`, the inverse friction constant for each group is
98       set with :mdp:`tau-t`.  The parameters :mdp:`tcoupl` and :mdp:`nsttcouple`
99       are ignored. The random generator is initialized with
100       :mdp:`ld-seed`. When used as a thermostat, an appropriate value
101       for :mdp:`tau-t` is 2 ps, since this results in a friction that
102       is lower than the internal friction of water, while it is high
103       enough to remove excess heat NOTE: temperature deviations decay
104       twice as fast as with a Berendsen thermostat with the same
105       :mdp:`tau-t`.
107    .. mdp-value:: bd
109       An Euler integrator for Brownian or position Langevin dynamics,
110       the velocity is the force divided by a friction coefficient
111       (:mdp:`bd-fric`) plus random thermal noise (:mdp:`ref-t`). When
112       :mdp:`bd-fric` is 0, the friction coefficient for each particle
113       is calculated as mass/ :mdp:`tau-t`, as for the integrator
114       :mdp-value:`integrator=sd`. The random generator is initialized
115       with :mdp:`ld-seed`.
117    .. mdp-value:: steep
119       A steepest descent algorithm for energy minimization. The
120       maximum step size is :mdp:`emstep`, the tolerance is
121       :mdp:`emtol`.
123    .. mdp-value:: cg
125       A conjugate gradient algorithm for energy minimization, the
126       tolerance is :mdp:`emtol`. CG is more efficient when a steepest
127       descent step is done every once in a while, this is determined
128       by :mdp:`nstcgsteep`. For a minimization prior to a normal mode
129       analysis, which requires a very high accuracy, |Gromacs| should be
130       compiled in double precision.
132    .. mdp-value:: l-bfgs
134       A quasi-Newtonian algorithm for energy minimization according to
135       the low-memory Broyden-Fletcher-Goldfarb-Shanno approach. In
136       practice this seems to converge faster than Conjugate Gradients,
137       but due to the correction steps necessary it is not (yet)
138       parallelized.
140    .. mdp-value:: nm
142       Normal mode analysis is performed on the structure in the :ref:`tpr`
143       file.  |Gromacs| should be compiled in double precision.
145    .. mdp-value:: tpi
147       Test particle insertion. The last molecule in the topology is
148       the test particle. A trajectory must be provided to ``mdrun
149       -rerun``. This trajectory should not contain the molecule to be
150       inserted. Insertions are performed :mdp:`nsteps` times in each
151       frame at random locations and with random orientiations of the
152       molecule. When :mdp:`nstlist` is larger than one,
153       :mdp:`nstlist` insertions are performed in a sphere with radius
154       :mdp:`rtpi` around a the same random location using the same
155       pair list. Since pair list construction is expensive,
156       one can perform several extra insertions with the same list
157       almost for free. The random seed is set with
158       :mdp:`ld-seed`. The temperature for the Boltzmann weighting is
159       set with :mdp:`ref-t`, this should match the temperature of the
160       simulation of the original trajectory. Dispersion correction is
161       implemented correctly for TPI. All relevant quantities are
162       written to the file specified with ``mdrun -tpi``. The
163       distribution of insertion energies is written to the file
164       specified with ``mdrun -tpid``. No trajectory or energy file is
165       written. Parallel TPI gives identical results to single-node
166       TPI. For charged molecules, using PME with a fine grid is most
167       accurate and also efficient, since the potential in the system
168       only needs to be calculated once per frame.
170    .. mdp-value:: tpic
172       Test particle insertion into a predefined cavity location. The
173       procedure is the same as for :mdp-value:`integrator=tpi`, except
174       that one coordinate extra is read from the trajectory, which is
175       used as the insertion location. The molecule to be inserted
176       should be centered at 0,0,0. |Gromacs| does not do this for you,
177       since for different situations a different way of centering
178       might be optimal. Also :mdp:`rtpi` sets the radius for the
179       sphere around this location. Neighbor searching is done only
180       once per frame, :mdp:`nstlist` is not used. Parallel
181       :mdp-value:`integrator=tpic` gives identical results to
182       single-rank :mdp-value:`integrator=tpic`.
184    .. mdp-value:: mimic
186       Enable MiMiC QM/MM coupling to run hybrid molecular dynamics.
187       Keey in mind that its required to launch CPMD compiled with MiMiC as well.
188       In this mode all options regarding integration (T-coupling, P-coupling,
189       timestep and number of steps) are ignored as CPMD will do the integration
190       instead. Options related to forces computation (cutoffs, PME parameters,
191       etc.) are working as usual. Atom selection to define QM atoms is read
192       from :mdp:`QMMM-grps`
194 .. mdp:: tinit
196         (0) [ps]
197         starting time for your run (only makes sense for time-based
198         integrators)
200 .. mdp:: dt
202         (0.001) [ps]
203         time step for integration (only makes sense for time-based
204         integrators)
206 .. mdp:: nsteps
208         (0)
209         maximum number of steps to integrate or minimize, -1 is no
210         maximum
212 .. mdp:: init-step
214         (0)
215         The starting step. The time at step i in a run is
216         calculated as: t = :mdp:`tinit` + :mdp:`dt` *
217         (:mdp:`init-step` + i). The free-energy lambda is calculated
218         as: lambda = :mdp:`init-lambda` + :mdp:`delta-lambda` *
219         (:mdp:`init-step` + i). Also non-equilibrium MD parameters can
220         depend on the step number. Thus for exact restarts or redoing
221         part of a run it might be necessary to set :mdp:`init-step` to
222         the step number of the restart frame. :ref:`gmx convert-tpr`
223         does this automatically.
225 .. mdp:: simulation-part
227          (0)
228          A simulation can consist of multiple parts, each of which has
229          a part number. This option specifies what that number will
230          be, which helps keep track of parts that are logically the
231          same simulation. This option is generally useful to set only
232          when coping with a crashed simulation where files were lost.
234 .. mdp:: comm-mode
236    .. mdp-value:: Linear
238       Remove center of mass translational velocity
240    .. mdp-value:: Angular
242       Remove center of mass translational and rotational velocity
244    .. mdp-value:: Linear-acceleration-correction
246       Remove center of mass translational velocity. Correct the center of
247       mass position assuming linear acceleration over :mdp:`nstcomm` steps.
248       This is useful for cases where an acceleration is expected on the
249       center of mass which is nearly constant over :mdp:`nstcomm` steps.
250       This can occur for example when pulling on a group using an absolute
251       reference.
253    .. mdp-value:: None
255       No restriction on the center of mass motion
257 .. mdp:: nstcomm
259    (100) [steps]
260    frequency for center of mass motion removal
262 .. mdp:: comm-grps
264    group(s) for center of mass motion removal, default is the whole
265    system
268 Langevin dynamics
269 ^^^^^^^^^^^^^^^^^
271 .. mdp:: bd-fric
273    (0) [amu ps\ :sup:`-1`]
274    Brownian dynamics friction coefficient. When :mdp:`bd-fric` is 0,
275    the friction coefficient for each particle is calculated as mass/
276    :mdp:`tau-t`.
278 .. mdp:: ld-seed
280    (-1) [integer]
281    used to initialize random generator for thermal noise for
282    stochastic and Brownian dynamics. When :mdp:`ld-seed` is set to -1,
283    a pseudo random seed is used. When running BD or SD on multiple
284    processors, each processor uses a seed equal to :mdp:`ld-seed` plus
285    the processor number.
288 Energy minimization
289 ^^^^^^^^^^^^^^^^^^^
291 .. mdp:: emtol
293    (10.0) [kJ mol\ :sup:`-1` nm\ :sup:`-1`]
294    the minimization is converged when the maximum force is smaller
295    than this value
297 .. mdp:: emstep
299    (0.01) [nm]
300    initial step-size
302 .. mdp:: nstcgsteep
304    (1000) [steps]
305    frequency of performing 1 steepest descent step while doing
306    conjugate gradient energy minimization.
308 .. mdp:: nbfgscorr
310    (10)
311    Number of correction steps to use for L-BFGS minimization. A higher
312    number is (at least theoretically) more accurate, but slower.
315 Shell Molecular Dynamics
316 ^^^^^^^^^^^^^^^^^^^^^^^^
318 When shells or flexible constraints are present in the system the
319 positions of the shells and the lengths of the flexible constraints
320 are optimized at every time step until either the RMS force on the
321 shells and constraints is less than :mdp:`emtol`, or a maximum number
322 of iterations :mdp:`niter` has been reached. Minimization is converged
323 when the maximum force is smaller than :mdp:`emtol`. For shell MD this
324 value should be 1.0 at most.
326 .. mdp:: niter
328    (20)
329    maximum number of iterations for optimizing the shell positions and
330    the flexible constraints.
332 .. mdp:: fcstep
334    (0) [ps\ :sup:`2`]
335    the step size for optimizing the flexible constraints. Should be
336    chosen as mu/(d2V/dq2) where mu is the reduced mass of two
337    particles in a flexible constraint and d2V/dq2 is the second
338    derivative of the potential in the constraint direction. Hopefully
339    this number does not differ too much between the flexible
340    constraints, as the number of iterations and thus the runtime is
341    very sensitive to fcstep. Try several values!
344 Test particle insertion
345 ^^^^^^^^^^^^^^^^^^^^^^^
347 .. mdp:: rtpi
349    (0.05) [nm]
350    the test particle insertion radius, see integrators
351    :mdp-value:`integrator=tpi` and :mdp-value:`integrator=tpic`
354 Output control
355 ^^^^^^^^^^^^^^
357 .. mdp:: nstxout
359    (0) [steps]
360    number of steps that elapse between writing coordinates to the output
361    trajectory file (:ref:`trr`), the last coordinates are always written
362    unless 0, which means coordinates are not written into the trajectory
363    file.
365 .. mdp:: nstvout
367    (0) [steps]
368    number of steps that elapse between writing velocities to the output
369    trajectory file (:ref:`trr`), the last velocities are always written
370    unless 0, which means velocities are not written into the trajectory
371    file.
373 .. mdp:: nstfout
375    (0) [steps]
376    number of steps that elapse between writing forces to the output
377    trajectory file (:ref:`trr`), the last forces are always written,
378    unless 0, which means forces are not written into the trajectory
379    file.
381 .. mdp:: nstlog
383    (1000) [steps]
384    number of steps that elapse between writing energies to the log
385    file, the last energies are always written.
387 .. mdp:: nstcalcenergy
389    (100)
390    number of steps that elapse between calculating the energies, 0 is
391    never. This option is only relevant with dynamics. This option affects the
392    performance in parallel simulations, because calculating energies
393    requires global communication between all processes which can
394    become a bottleneck at high parallelization.
396 .. mdp:: nstenergy
398    (1000) [steps]
399    number of steps that elapse between writing energies to energy file,
400    the last energies are always written, should be a multiple of
401    :mdp:`nstcalcenergy`. Note that the exact sums and fluctuations
402    over all MD steps modulo :mdp:`nstcalcenergy` are stored in the
403    energy file, so :ref:`gmx energy` can report exact energy averages
404    and fluctuations also when :mdp:`nstenergy` > 1
406 .. mdp:: nstxout-compressed
408    (0) [steps]
409    number of steps that elapse between writing position coordinates
410    using lossy compression (:ref:`xtc` file), 0 for not writing
411    compressed coordinates output.
413 .. mdp:: compressed-x-precision
415    (1000) [real]
416    precision with which to write to the compressed trajectory file
418 .. mdp:: compressed-x-grps
420    group(s) to write to the compressed trajectory file, by default the
421    whole system is written (if :mdp:`nstxout-compressed` > 0)
423 .. mdp:: energygrps
425    group(s) for which to write to write short-ranged non-bonded
426    potential energies to the energy file (not supported on GPUs)
429 Neighbor searching
430 ^^^^^^^^^^^^^^^^^^
432 .. mdp:: cutoff-scheme
434    .. mdp-value:: Verlet
436       Generate a pair list with buffering. The buffer size is
437       automatically set based on :mdp:`verlet-buffer-tolerance`,
438       unless this is set to -1, in which case :mdp:`rlist` will be
439       used.
441    .. mdp-value:: group
443       Generate a pair list for groups of atoms, corresponding
444       to the charge groups in the topology. This option is no longer
445       supported.
447 .. mdp:: nstlist
449    (10) [steps]
451    .. mdp-value:: >0
453       Frequency to update the neighbor list. When dynamics and
454       :mdp:`verlet-buffer-tolerance` set, :mdp:`nstlist` is actually
455       a minimum value and :ref:`gmx mdrun` might increase it, unless
456       it is set to 1. With parallel simulations and/or non-bonded
457       force calculation on the GPU, a value of 20 or 40 often gives
458       the best performance.
460    .. mdp-value:: 0
462       The neighbor list is only constructed once and never
463       updated. This is mainly useful for vacuum simulations in which
464       all particles see each other. But vacuum simulations are
465       (temporarily) not supported.
467    .. mdp-value:: <0
469       Unused.
471 .. mdp:: pbc
473    .. mdp-value:: xyz
475       Use periodic boundary conditions in all directions.
477    .. mdp-value:: no
479       Use no periodic boundary conditions, ignore the box. To simulate
480       without cut-offs, set all cut-offs and :mdp:`nstlist` to 0. For
481       best performance without cut-offs on a single MPI rank, set
482       :mdp:`nstlist` to zero and :mdp-value:`ns-type=simple`.
484    .. mdp-value:: xy
486       Use periodic boundary conditions in x and y directions
487       only. This works only with :mdp-value:`ns-type=grid` and can be used
488       in combination with walls_. Without walls or with only one wall
489       the system size is infinite in the z direction. Therefore
490       pressure coupling or Ewald summation methods can not be
491       used. These disadvantages do not apply when two walls are used.
493 .. mdp:: periodic-molecules
495    .. mdp-value:: no
497       molecules are finite, fast molecular PBC can be used
499    .. mdp-value:: yes
501       for systems with molecules that couple to themselves through the
502       periodic boundary conditions, this requires a slower PBC
503       algorithm and molecules are not made whole in the output
505 .. mdp:: verlet-buffer-tolerance
507    (0.005) [kJ mol\ :sup:`-1` ps\ :sup:`-1`]
509    Used when performing a simulation with dynamics. This sets
510    the maximum allowed error for pair interactions per particle caused
511    by the Verlet buffer, which indirectly sets :mdp:`rlist`. As both
512    :mdp:`nstlist` and the Verlet buffer size are fixed (for
513    performance reasons), particle pairs not in the pair list can
514    occasionally get within the cut-off distance during
515    :mdp:`nstlist` -1 steps. This causes very small jumps in the
516    energy. In a constant-temperature ensemble, these very small energy
517    jumps can be estimated for a given cut-off and :mdp:`rlist`. The
518    estimate assumes a homogeneous particle distribution, hence the
519    errors might be slightly underestimated for multi-phase
520    systems. (See the `reference manual`_ for details). For longer
521    pair-list life-time (:mdp:`nstlist` -1) * :mdp:`dt` the buffer is
522    overestimated, because the interactions between particles are
523    ignored. Combined with cancellation of errors, the actual drift of
524    the total energy is usually one to two orders of magnitude
525    smaller. Note that the generated buffer size takes into account
526    that the |Gromacs| pair-list setup leads to a reduction in the
527    drift by a factor 10, compared to a simple particle-pair based
528    list. Without dynamics (energy minimization etc.), the buffer is 5%
529    of the cut-off. For NVE simulations the initial temperature is
530    used, unless this is zero, in which case a buffer of 10% is
531    used. For NVE simulations the tolerance usually needs to be lowered
532    to achieve proper energy conservation on the nanosecond time
533    scale. To override the automated buffer setting, use
534    :mdp:`verlet-buffer-tolerance` =-1 and set :mdp:`rlist` manually.
536 .. mdp:: rlist
538    (1) [nm]
539    Cut-off distance for the short-range neighbor list. With dynamics,
540    this is by default set by the :mdp:`verlet-buffer-tolerance` option
541    and the value of :mdp:`rlist` is ignored. Without dynamics, this
542    is by default set to the maximum cut-off plus 5% buffer, except
543    for test particle insertion, where the buffer is managed exactly
544    and automatically. For NVE simulations, where the automated
545    setting is not possible, the advised procedure is to run :ref:`gmx grompp`
546    with an NVT setup with the expected temperature and copy the resulting
547    value of :mdp:`rlist` to the NVE setup.
550 Electrostatics
551 ^^^^^^^^^^^^^^
553 .. mdp:: coulombtype
555    .. mdp-value:: Cut-off
557       Plain cut-off with pair list radius :mdp:`rlist` and
558       Coulomb cut-off :mdp:`rcoulomb`, where :mdp:`rlist` >=
559       :mdp:`rcoulomb`.
561    .. mdp-value:: Ewald
563       Classical Ewald sum electrostatics. The real-space cut-off
564       :mdp:`rcoulomb` should be equal to :mdp:`rlist`. Use *e.g.*
565       :mdp:`rlist` =0.9, :mdp:`rcoulomb` =0.9. The highest magnitude
566       of wave vectors used in reciprocal space is controlled by
567       :mdp:`fourierspacing`. The relative accuracy of
568       direct/reciprocal space is controlled by :mdp:`ewald-rtol`.
570       NOTE: Ewald scales as O(N\ :sup:`3/2`) and is thus extremely slow for
571       large systems. It is included mainly for reference - in most
572       cases PME will perform much better.
574    .. mdp-value:: PME
576       Fast smooth Particle-Mesh Ewald (SPME) electrostatics. Direct
577       space is similar to the Ewald sum, while the reciprocal part is
578       performed with FFTs. Grid dimensions are controlled with
579       :mdp:`fourierspacing` and the interpolation order with
580       :mdp:`pme-order`. With a grid spacing of 0.1 nm and cubic
581       interpolation the electrostatic forces have an accuracy of
582       2-3*10\ :sup:`-4`. Since the error from the vdw-cutoff is larger than
583       this you might try 0.15 nm. When running in parallel the
584       interpolation parallelizes better than the FFT, so try
585       decreasing grid dimensions while increasing interpolation.
587    .. mdp-value:: P3M-AD
589       Particle-Particle Particle-Mesh algorithm with analytical
590       derivative for for long range electrostatic interactions. The
591       method and code is identical to SPME, except that the influence
592       function is optimized for the grid. This gives a slight increase
593       in accuracy.
595    .. mdp-value:: Reaction-Field
597       Reaction field electrostatics with Coulomb cut-off
598       :mdp:`rcoulomb`, where :mdp:`rlist` >= :mdp:`rvdw`. The
599       dielectric constant beyond the cut-off is
600       :mdp:`epsilon-rf`. The dielectric constant can be set to
601       infinity by setting :mdp:`epsilon-rf` =0.
603    .. mdp-value:: User
605       Currently unsupported.
606       :ref:`gmx mdrun` will now expect to find a file ``table.xvg``
607       with user-defined potential functions for repulsion, dispersion
608       and Coulomb. When pair interactions are present, :ref:`gmx
609       mdrun` also expects to find a file ``tablep.xvg`` for the pair
610       interactions. When the same interactions should be used for
611       non-bonded and pair interactions the user can specify the same
612       file name for both table files. These files should contain 7
613       columns: the ``x`` value, ``f(x)``, ``-f'(x)``, ``g(x)``,
614       ``-g'(x)``, ``h(x)``, ``-h'(x)``, where ``f(x)`` is the Coulomb
615       function, ``g(x)`` the dispersion function and ``h(x)`` the
616       repulsion function. When :mdp:`vdwtype` is not set to User the
617       values for ``g``, ``-g'``, ``h`` and ``-h'`` are ignored. For
618       the non-bonded interactions ``x`` values should run from 0 to
619       the largest cut-off distance + :mdp:`table-extension` and
620       should be uniformly spaced. For the pair interactions the table
621       length in the file will be used. The optimal spacing, which is
622       used for non-user tables, is ``0.002 nm`` when you run in mixed
623       precision or ``0.0005 nm`` when you run in double precision. The
624       function value at ``x=0`` is not important. More information is
625       in the printed manual.
627    .. mdp-value:: PME-Switch
629       Currently unsupported.
630       A combination of PME and a switch function for the direct-space
631       part (see above). :mdp:`rcoulomb` is allowed to be smaller than
632       :mdp:`rlist`.
634    .. mdp-value:: PME-User
636       Currently unsupported.
637       A combination of PME and user tables (see
638       above). :mdp:`rcoulomb` is allowed to be smaller than
639       :mdp:`rlist`. The PME mesh contribution is subtracted from the
640       user table by :ref:`gmx mdrun`. Because of this subtraction the
641       user tables should contain about 10 decimal places.
643    .. mdp-value:: PME-User-Switch
645       Currently unsupported.
646       A combination of PME-User and a switching function (see
647       above). The switching function is applied to final
648       particle-particle interaction, *i.e.* both to the user supplied
649       function and the PME Mesh correction part.
651 .. mdp:: coulomb-modifier
653    .. mdp-value:: Potential-shift
655       Shift the Coulomb potential by a constant such that it is zero
656       at the cut-off. This makes the potential the integral of the
657       force. Note that this does not affect the forces or the
658       sampling.
660    .. mdp-value:: None
662       Use an unmodified Coulomb potential. This can be useful
663       when comparing energies with those computed with other software.
665 .. mdp:: rcoulomb-switch
667    (0) [nm]
668    where to start switching the Coulomb potential, only relevant
669    when force or potential switching is used
671 .. mdp:: rcoulomb
673    (1) [nm]
674    The distance for the Coulomb cut-off. Note that with PME this value
675    can be increased by the PME tuning in :ref:`gmx mdrun` along with
676    the PME grid spacing.
678 .. mdp:: epsilon-r
680    (1)
681    The relative dielectric constant. A value of 0 means infinity.
683 .. mdp:: epsilon-rf
685    (0)
686    The relative dielectric constant of the reaction field. This
687    is only used with reaction-field electrostatics. A value of 0
688    means infinity.
691 Van der Waals
692 ^^^^^^^^^^^^^
694 .. mdp:: vdwtype
696    .. mdp-value:: Cut-off
698       Plain cut-off with pair list radius :mdp:`rlist` and VdW
699       cut-off :mdp:`rvdw`, where :mdp:`rlist` >= :mdp:`rvdw`.
701    .. mdp-value:: PME
703       Fast smooth Particle-mesh Ewald (SPME) for VdW interactions. The
704       grid dimensions are controlled with :mdp:`fourierspacing` in
705       the same way as for electrostatics, and the interpolation order
706       is controlled with :mdp:`pme-order`. The relative accuracy of
707       direct/reciprocal space is controlled by :mdp:`ewald-rtol-lj`,
708       and the specific combination rules that are to be used by the
709       reciprocal routine are set using :mdp:`lj-pme-comb-rule`.
711    .. mdp-value:: Shift
713       This functionality is deprecated and replaced by using
714       :mdp-value:`vdwtype=Cut-off` with :mdp-value:`vdw-modifier=Force-switch`.
715       The LJ (not Buckingham) potential is decreased over the whole range and
716       the forces decay smoothly to zero between :mdp:`rvdw-switch` and
717       :mdp:`rvdw`.
719    .. mdp-value:: Switch
721       This functionality is deprecated and replaced by using
722       :mdp-value:`vdwtype=Cut-off` with :mdp-value:`vdw-modifier=Potential-switch`.
723       The LJ (not Buckingham) potential is normal out to :mdp:`rvdw-switch`, after
724       which it is switched off to reach zero at :mdp:`rvdw`. Both the
725       potential and force functions are continuously smooth, but be
726       aware that all switch functions will give rise to a bulge
727       (increase) in the force (since we are switching the
728       potential).
730    .. mdp-value:: User
732       Currently unsupported.
733       See user for :mdp:`coulombtype`. The function value at zero is
734       not important. When you want to use LJ correction, make sure
735       that :mdp:`rvdw` corresponds to the cut-off in the user-defined
736       function. When :mdp:`coulombtype` is not set to User the values
737       for the ``f`` and ``-f'`` columns are ignored.
739 .. mdp:: vdw-modifier
741    .. mdp-value:: Potential-shift
743       Shift the Van der Waals potential by a constant such that it is
744       zero at the cut-off. This makes the potential the integral of
745       the force. Note that this does not affect the forces or the
746       sampling.
748    .. mdp-value:: None
750       Use an unmodified Van der Waals potential. This can be useful
751       when comparing energies with those computed with other software.
753    .. mdp-value:: Force-switch
755       Smoothly switches the forces to zero between :mdp:`rvdw-switch`
756       and :mdp:`rvdw`. This shifts the potential shift over the whole
757       range and switches it to zero at the cut-off. Note that this is
758       more expensive to calculate than a plain cut-off and it is not
759       required for energy conservation, since Potential-shift
760       conserves energy just as well.
762    .. mdp-value:: Potential-switch
764       Smoothly switches the potential to zero between
765       :mdp:`rvdw-switch` and :mdp:`rvdw`. Note that this introduces
766       articifically large forces in the switching region and is much
767       more expensive to calculate. This option should only be used if
768       the force field you are using requires this.
770 .. mdp:: rvdw-switch
772    (0) [nm]
773    where to start switching the LJ force and possibly the potential,
774    only relevant when force or potential switching is used
776 .. mdp:: rvdw
778    (1) [nm]
779    distance for the LJ or Buckingham cut-off
781 .. mdp:: DispCorr
783    .. mdp-value:: no
785       don't apply any correction
787    .. mdp-value:: EnerPres
789       apply long range dispersion corrections for Energy and Pressure
791    .. mdp-value:: Ener
793       apply long range dispersion corrections for Energy only
796 Tables
797 ^^^^^^
799 .. mdp:: table-extension
801    (1) [nm]
802    Extension of the non-bonded potential lookup tables beyond the
803    largest cut-off distance. With actual non-bonded interactions
804    the tables are never accessed beyond the cut-off. But a longer
805    table length might be needed for the 1-4 interactions, which
806    are always tabulated irrespective of the use of tables for
807    the non-bonded interactions.
809 .. mdp:: energygrp-table
811    Currently unsupported.
812    When user tables are used for electrostatics and/or VdW, here one
813    can give pairs of energy groups for which seperate user tables
814    should be used. The two energy groups will be appended to the table
815    file name, in order of their definition in :mdp:`energygrps`,
816    seperated by underscores. For example, if ``energygrps = Na Cl
817    Sol`` and ``energygrp-table = Na Na Na Cl``, :ref:`gmx mdrun` will
818    read ``table_Na_Na.xvg`` and ``table_Na_Cl.xvg`` in addition to the
819    normal ``table.xvg`` which will be used for all other energy group
820    pairs.
823 Ewald
824 ^^^^^
826 .. mdp:: fourierspacing
828    (0.12) [nm]
829    For ordinary Ewald, the ratio of the box dimensions and the spacing
830    determines a lower bound for the number of wave vectors to use in
831    each (signed) direction. For PME and P3M, that ratio determines a
832    lower bound for the number of Fourier-space grid points that will
833    be used along that axis. In all cases, the number for each
834    direction can be overridden by entering a non-zero value for that
835    :mdp:`fourier-nx` direction. For optimizing the relative load of
836    the particle-particle interactions and the mesh part of PME, it is
837    useful to know that the accuracy of the electrostatics remains
838    nearly constant when the Coulomb cut-off and the PME grid spacing
839    are scaled by the same factor. Note that this spacing can be scaled
840    up along with :mdp:`rcoulomb` by the PME tuning in :ref:`gmx mdrun`.
842 .. mdp:: fourier-nx
843 .. mdp:: fourier-ny
844 .. mdp:: fourier-nz
846    (0)
847    Highest magnitude of wave vectors in reciprocal space when using Ewald.
848    Grid size when using PME or P3M. These values override
849    :mdp:`fourierspacing` per direction. The best choice is powers of
850    2, 3, 5 and 7. Avoid large primes. Note that these grid sizes can
851    be reduced along with scaling up :mdp:`rcoulomb` by the PME tuning
852    in :ref:`gmx mdrun`.
854 .. mdp:: pme-order
856    (4)
857    Interpolation order for PME. 4 equals cubic interpolation. You
858    might try 6/8/10 when running in parallel and simultaneously
859    decrease grid dimension.
861 .. mdp:: ewald-rtol
863    (10\ :sup:`-5`)
864    The relative strength of the Ewald-shifted direct potential at
865    :mdp:`rcoulomb` is given by :mdp:`ewald-rtol`. Decreasing this
866    will give a more accurate direct sum, but then you need more wave
867    vectors for the reciprocal sum.
869 .. mdp:: ewald-rtol-lj
871    (10\ :sup:`-3`)
872    When doing PME for VdW-interactions, :mdp:`ewald-rtol-lj` is used
873    to control the relative strength of the dispersion potential at
874    :mdp:`rvdw` in the same way as :mdp:`ewald-rtol` controls the
875    electrostatic potential.
877 .. mdp:: lj-pme-comb-rule
879    (Geometric)
880    The combination rules used to combine VdW-parameters in the
881    reciprocal part of LJ-PME. Geometric rules are much faster than
882    Lorentz-Berthelot and usually the recommended choice, even when the
883    rest of the force field uses the Lorentz-Berthelot rules.
885    .. mdp-value:: Geometric
887       Apply geometric combination rules
889    .. mdp-value:: Lorentz-Berthelot
891       Apply Lorentz-Berthelot combination rules
893 .. mdp:: ewald-geometry
895    .. mdp-value:: 3d
897       The Ewald sum is performed in all three dimensions.
899    .. mdp-value:: 3dc
901       The reciprocal sum is still performed in 3D, but a force and
902       potential correction applied in the ``z`` dimension to produce a
903       pseudo-2D summation. If your system has a slab geometry in the
904       ``x-y`` plane you can try to increase the ``z``-dimension of the box
905       (a box height of 3 times the slab height is usually ok) and use
906       this option.
908 .. mdp:: epsilon-surface
910    (0)
911    This controls the dipole correction to the Ewald summation in
912    3D. The default value of zero means it is turned off. Turn it on by
913    setting it to the value of the relative permittivity of the
914    imaginary surface around your infinite system. Be careful - you
915    shouldn't use this if you have free mobile charges in your
916    system. This value does not affect the slab 3DC variant of the long
917    range corrections.
920 Temperature coupling
921 ^^^^^^^^^^^^^^^^^^^^
923 .. mdp:: tcoupl
925    .. mdp-value:: no
927       No temperature coupling.
929    .. mdp-value:: berendsen
931       Temperature coupling with a Berendsen thermostat to a bath with
932       temperature :mdp:`ref-t`, with time constant
933       :mdp:`tau-t`. Several groups can be coupled separately, these
934       are specified in the :mdp:`tc-grps` field separated by spaces.
936    .. mdp-value:: nose-hoover
938       Temperature coupling using a Nose-Hoover extended ensemble. The
939       reference temperature and coupling groups are selected as above,
940       but in this case :mdp:`tau-t` controls the period of the
941       temperature fluctuations at equilibrium, which is slightly
942       different from a relaxation time. For NVT simulations the
943       conserved energy quantity is written to the energy and log files.
945    .. mdp-value:: andersen
947       Temperature coupling by randomizing a fraction of the particle velocities
948       at each timestep. Reference temperature and coupling groups are
949       selected as above. :mdp:`tau-t` is the average time between
950       randomization of each molecule. Inhibits particle dynamics
951       somewhat, but little or no ergodicity issues. Currently only
952       implemented with velocity Verlet, and not implemented with
953       constraints.
955    .. mdp-value:: andersen-massive
957       Temperature coupling by randomizing velocities of all particles at
958       infrequent timesteps. Reference temperature and coupling groups are
959       selected as above. :mdp:`tau-t` is the time between
960       randomization of all molecules. Inhibits particle dynamics
961       somewhat, but little or no ergodicity issues. Currently only
962       implemented with velocity Verlet.
964    .. mdp-value:: v-rescale
966       Temperature coupling using velocity rescaling with a stochastic
967       term (JCP 126, 014101). This thermostat is similar to Berendsen
968       coupling, with the same scaling using :mdp:`tau-t`, but the
969       stochastic term ensures that a proper canonical ensemble is
970       generated. The random seed is set with :mdp:`ld-seed`. This
971       thermostat works correctly even for :mdp:`tau-t` =0. For NVT
972       simulations the conserved energy quantity is written to the
973       energy and log file.
975 .. mdp:: nsttcouple
977    (-1)
978    The frequency for coupling the temperature. The default value of -1
979    sets :mdp:`nsttcouple` equal to :mdp:`nstlist`, unless
980    :mdp:`nstlist` <=0, then a value of 10 is used. For velocity
981    Verlet integrators :mdp:`nsttcouple` is set to 1.
983 .. mdp:: nh-chain-length
985    (10)
986    The number of chained Nose-Hoover thermostats for velocity Verlet
987    integrators, the leap-frog :mdp-value:`integrator=md` integrator
988    only supports 1. Data for the NH chain variables is not printed
989    to the :ref:`edr` file by default, but can be turned on with the
990    :mdp:`print-nose-hoover-chain-variables` option.
992 .. mdp:: print-nose-hoover-chain-variables
994    .. mdp-value:: no
996       Do not store Nose-Hoover chain variables in the energy file.
998    .. mdp-value:: yes
1000       Store all positions and velocities of the Nose-Hoover chain
1001       in the energy file.
1003 .. mdp:: tc-grps
1005    groups to couple to separate temperature baths
1007 .. mdp:: tau-t
1009    [ps]
1010    time constant for coupling (one for each group in
1011    :mdp:`tc-grps`), -1 means no temperature coupling
1013 .. mdp:: ref-t
1015    [K]
1016    reference temperature for coupling (one for each group in
1017    :mdp:`tc-grps`)
1020 Pressure coupling
1021 ^^^^^^^^^^^^^^^^^
1023 .. mdp:: pcoupl
1025    .. mdp-value:: no
1027       No pressure coupling. This means a fixed box size.
1029    .. mdp-value:: Berendsen
1031       Exponential relaxation pressure coupling with time constant
1032       :mdp:`tau-p`. The box is scaled every :mdp:`nstpcouple` steps. It has been
1033       argued that this does not yield a correct thermodynamic
1034       ensemble, but it is the most efficient way to scale a box at the
1035       beginning of a run.
1037    .. mdp-value:: Parrinello-Rahman
1039       Extended-ensemble pressure coupling where the box vectors are
1040       subject to an equation of motion. The equation of motion for the
1041       atoms is coupled to this. No instantaneous scaling takes
1042       place. As for Nose-Hoover temperature coupling the time constant
1043       :mdp:`tau-p` is the period of pressure fluctuations at
1044       equilibrium. This is probably a better method when you want to
1045       apply pressure scaling during data collection, but beware that
1046       you can get very large oscillations if you are starting from a
1047       different pressure. For simulations where the exact fluctations
1048       of the NPT ensemble are important, or if the pressure coupling
1049       time is very short it may not be appropriate, as the previous
1050       time step pressure is used in some steps of the |Gromacs|
1051       implementation for the current time step pressure.
1053    .. mdp-value:: MTTK
1055       Martyna-Tuckerman-Tobias-Klein implementation, only useable with
1056       :mdp-value:`integrator=md-vv` or :mdp-value:`integrator=md-vv-avek`, very similar to
1057       Parrinello-Rahman. As for Nose-Hoover temperature coupling the
1058       time constant :mdp:`tau-p` is the period of pressure
1059       fluctuations at equilibrium. This is probably a better method
1060       when you want to apply pressure scaling during data collection,
1061       but beware that you can get very large oscillations if you are
1062       starting from a different pressure. Currently (as of version
1063       5.1), it only supports isotropic scaling, and only works without
1064       constraints.
1066 .. mdp:: pcoupltype
1068    Specifies the kind of isotropy of the pressure coupling used. Each
1069    kind takes one or more values for :mdp:`compressibility` and
1070    :mdp:`ref-p`. Only a single value is permitted for :mdp:`tau-p`.
1072    .. mdp-value:: isotropic
1074       Isotropic pressure coupling with time constant
1075       :mdp:`tau-p`. One value each for :mdp:`compressibility` and
1076       :mdp:`ref-p` is required.
1078    .. mdp-value:: semiisotropic
1080       Pressure coupling which is isotropic in the ``x`` and ``y``
1081       direction, but different in the ``z`` direction. This can be
1082       useful for membrane simulations. Two values each for
1083       :mdp:`compressibility` and :mdp:`ref-p` are required, for
1084       ``x/y`` and ``z`` directions respectively.
1086    .. mdp-value:: anisotropic
1088       Same as before, but 6 values are needed for ``xx``, ``yy``, ``zz``,
1089       ``xy/yx``, ``xz/zx`` and ``yz/zy`` components,
1090       respectively. When the off-diagonal compressibilities are set to
1091       zero, a rectangular box will stay rectangular. Beware that
1092       anisotropic scaling can lead to extreme deformation of the
1093       simulation box.
1095    .. mdp-value:: surface-tension
1097       Surface tension coupling for surfaces parallel to the
1098       xy-plane. Uses normal pressure coupling for the ``z``-direction,
1099       while the surface tension is coupled to the ``x/y`` dimensions of
1100       the box. The first :mdp:`ref-p` value is the reference surface
1101       tension times the number of surfaces ``bar nm``, the second
1102       value is the reference ``z``-pressure ``bar``. The two
1103       :mdp:`compressibility` values are the compressibility in the
1104       ``x/y`` and ``z`` direction respectively. The value for the
1105       ``z``-compressibility should be reasonably accurate since it
1106       influences the convergence of the surface-tension, it can also
1107       be set to zero to have a box with constant height.
1109 .. mdp:: nstpcouple
1111    (-1)
1112    The frequency for coupling the pressure. The default value of -1
1113    sets :mdp:`nstpcouple` equal to :mdp:`nstlist`, unless
1114    :mdp:`nstlist` <=0, then a value of 10 is used. For velocity
1115    Verlet integrators :mdp:`nstpcouple` is set to 1.
1117 .. mdp:: tau-p
1119    (1) [ps]
1120    The time constant for pressure coupling (one value for all
1121    directions).
1123 .. mdp:: compressibility
1125    [bar\ :sup:`-1`]
1126    The compressibility (NOTE: this is now really in bar\ :sup:`-1`) For water at 1
1127    atm and 300 K the compressibility is 4.5e-5 bar\ :sup:`-1`. The number of
1128    required values is implied by :mdp:`pcoupltype`.
1130 .. mdp:: ref-p
1132    [bar]
1133    The reference pressure for coupling. The number of required values
1134    is implied by :mdp:`pcoupltype`.
1136 .. mdp:: refcoord-scaling
1138    .. mdp-value:: no
1140       The reference coordinates for position restraints are not
1141       modified. Note that with this option the virial and pressure
1142       might be ill defined, see :ref:`here <reference-manual-position-restraints>`
1143       for more details.
1145    .. mdp-value:: all
1147       The reference coordinates are scaled with the scaling matrix of
1148       the pressure coupling.
1150    .. mdp-value:: com
1152       Scale the center of mass of the reference coordinates with the
1153       scaling matrix of the pressure coupling. The vectors of each
1154       reference coordinate to the center of mass are not scaled. Only
1155       one COM is used, even when there are multiple molecules with
1156       position restraints. For calculating the COM of the reference
1157       coordinates in the starting configuration, periodic boundary
1158       conditions are not taken into account. Note that with this option
1159       the virial and pressure might be ill defined, see
1160       :ref:`here <reference-manual-position-restraints>` for more details.
1163 Simulated annealing
1164 ^^^^^^^^^^^^^^^^^^^
1166 Simulated annealing is controlled separately for each temperature
1167 group in |Gromacs|. The reference temperature is a piecewise linear
1168 function, but you can use an arbitrary number of points for each
1169 group, and choose either a single sequence or a periodic behaviour for
1170 each group. The actual annealing is performed by dynamically changing
1171 the reference temperature used in the thermostat algorithm selected,
1172 so remember that the system will usually not instantaneously reach the
1173 reference temperature!
1175 .. mdp:: annealing
1177    Type of annealing for each temperature group
1179    .. mdp-value:: no
1181        No simulated annealing - just couple to reference temperature value.
1183    .. mdp-value:: single
1185        A single sequence of annealing points. If your simulation is
1186        longer than the time of the last point, the temperature will be
1187        coupled to this constant value after the annealing sequence has
1188        reached the last time point.
1190    .. mdp-value:: periodic
1192        The annealing will start over at the first reference point once
1193        the last reference time is reached. This is repeated until the
1194        simulation ends.
1196 .. mdp:: annealing-npoints
1198    A list with the number of annealing reference/control points used
1199    for each temperature group. Use 0 for groups that are not
1200    annealed. The number of entries should equal the number of
1201    temperature groups.
1203 .. mdp:: annealing-time
1205    List of times at the annealing reference/control points for each
1206    group. If you are using periodic annealing, the times will be used
1207    modulo the last value, *i.e.* if the values are 0, 5, 10, and 15,
1208    the coupling will restart at the 0ps value after 15ps, 30ps, 45ps,
1209    etc. The number of entries should equal the sum of the numbers
1210    given in :mdp:`annealing-npoints`.
1212 .. mdp:: annealing-temp
1214    List of temperatures at the annealing reference/control points for
1215    each group. The number of entries should equal the sum of the
1216    numbers given in :mdp:`annealing-npoints`.
1218 Confused? OK, let's use an example. Assume you have two temperature
1219 groups, set the group selections to ``annealing = single periodic``,
1220 the number of points of each group to ``annealing-npoints = 3 4``, the
1221 times to ``annealing-time = 0 3 6 0 2 4 6`` and finally temperatures
1222 to ``annealing-temp = 298 280 270 298 320 320 298``. The first group
1223 will be coupled to 298K at 0ps, but the reference temperature will
1224 drop linearly to reach 280K at 3ps, and then linearly between 280K and
1225 270K from 3ps to 6ps. After this is stays constant, at 270K. The
1226 second group is coupled to 298K at 0ps, it increases linearly to 320K
1227 at 2ps, where it stays constant until 4ps. Between 4ps and 6ps it
1228 decreases to 298K, and then it starts over with the same pattern
1229 again, *i.e.* rising linearly from 298K to 320K between 6ps and
1230 8ps. Check the summary printed by :ref:`gmx grompp` if you are unsure!
1233 Velocity generation
1234 ^^^^^^^^^^^^^^^^^^^
1236 .. mdp:: gen-vel
1238    .. mdp-value:: no
1240         Do not generate velocities. The velocities are set to zero
1241         when there are no velocities in the input structure file.
1243    .. mdp-value:: yes
1245         Generate velocities in :ref:`gmx grompp` according to a
1246         Maxwell distribution at temperature :mdp:`gen-temp`, with
1247         random seed :mdp:`gen-seed`. This is only meaningful with
1248         :mdp-value:`integrator=md`.
1250 .. mdp:: gen-temp
1252    (300) [K]
1253    temperature for Maxwell distribution
1255 .. mdp:: gen-seed
1257    (-1) [integer]
1258    used to initialize random generator for random velocities,
1259    when :mdp:`gen-seed` is set to -1, a pseudo random seed is
1260    used.
1263 Bonds
1264 ^^^^^
1266 .. mdp:: constraints
1268    Controls which bonds in the topology will be converted to rigid
1269    holonomic constraints. Note that typical rigid water models do not
1270    have bonds, but rather a specialized ``[settles]`` directive, so
1271    are not affected by this keyword.
1273    .. mdp-value:: none
1275       No bonds converted to constraints.
1277    .. mdp-value:: h-bonds
1279       Convert the bonds with H-atoms to constraints.
1281    .. mdp-value:: all-bonds
1283       Convert all bonds to constraints.
1285    .. mdp-value:: h-angles
1287       Convert all bonds to constraints and convert the angles that
1288       involve H-atoms to bond-constraints.
1290    .. mdp-value:: all-angles
1292       Convert all bonds to constraints and all angles to bond-constraints.
1294 .. mdp:: constraint-algorithm
1296    Chooses which solver satisfies any non-SETTLE holonomic
1297    constraints.
1299    .. mdp-value:: LINCS
1301       LINear Constraint Solver. With domain decomposition the parallel
1302       version P-LINCS is used. The accuracy in set with
1303       :mdp:`lincs-order`, which sets the number of matrices in the
1304       expansion for the matrix inversion. After the matrix inversion
1305       correction the algorithm does an iterative correction to
1306       compensate for lengthening due to rotation. The number of such
1307       iterations can be controlled with :mdp:`lincs-iter`. The root
1308       mean square relative constraint deviation is printed to the log
1309       file every :mdp:`nstlog` steps. If a bond rotates more than
1310       :mdp:`lincs-warnangle` in one step, a warning will be printed
1311       both to the log file and to ``stderr``. LINCS should not be used
1312       with coupled angle constraints.
1314    .. mdp-value:: SHAKE
1316       SHAKE is slightly slower and less stable than LINCS, but does
1317       work with angle constraints. The relative tolerance is set with
1318       :mdp:`shake-tol`, 0.0001 is a good value for "normal" MD. SHAKE
1319       does not support constraints between atoms on different
1320       decomposition domains, so it can only be used with domain
1321       decomposition when so-called update-groups are used, which is
1322       usally the case when only bonds involving hydrogens are
1323       constrained. SHAKE can not be used with energy minimization.
1325 .. mdp:: continuation
1327    This option was formerly known as ``unconstrained-start``.
1329    .. mdp-value:: no
1331       apply constraints to the start configuration and reset shells
1333    .. mdp-value:: yes
1335       do not apply constraints to the start configuration and do not
1336       reset shells, useful for exact coninuation and reruns
1338 .. mdp:: shake-tol
1340    (0.0001)
1341    relative tolerance for SHAKE
1343 .. mdp:: lincs-order
1345    (4)
1346    Highest order in the expansion of the constraint coupling
1347    matrix. When constraints form triangles, an additional expansion of
1348    the same order is applied on top of the normal expansion only for
1349    the couplings within such triangles. For "normal" MD simulations an
1350    order of 4 usually suffices, 6 is needed for large time-steps with
1351    virtual sites or BD. For accurate energy minimization an order of 8
1352    or more might be required. With domain decomposition, the cell size
1353    is limited by the distance spanned by :mdp:`lincs-order` +1
1354    constraints. When one wants to scale further than this limit, one
1355    can decrease :mdp:`lincs-order` and increase :mdp:`lincs-iter`,
1356    since the accuracy does not deteriorate when (1+ :mdp:`lincs-iter`
1357    )* :mdp:`lincs-order` remains constant.
1359 .. mdp:: lincs-iter
1361    (1)
1362    Number of iterations to correct for rotational lengthening in
1363    LINCS. For normal runs a single step is sufficient, but for NVE
1364    runs where you want to conserve energy accurately or for accurate
1365    energy minimization you might want to increase it to 2.
1367 .. mdp:: lincs-warnangle
1369    (30) [deg]
1370    maximum angle that a bond can rotate before LINCS will complain
1372 .. mdp:: morse
1374    .. mdp-value:: no
1376       bonds are represented by a harmonic potential
1378    .. mdp-value:: yes
1380       bonds are represented by a Morse potential
1383 Energy group exclusions
1384 ^^^^^^^^^^^^^^^^^^^^^^^
1386 .. mdp:: energygrp-excl
1388    Pairs of energy groups for which all non-bonded interactions are
1389    excluded. An example: if you have two energy groups ``Protein`` and
1390    ``SOL``, specifying ``energygrp-excl = Protein Protein SOL SOL``
1391    would give only the non-bonded interactions between the protein and
1392    the solvent. This is especially useful for speeding up energy
1393    calculations with ``mdrun -rerun`` and for excluding interactions
1394    within frozen groups.
1397 Walls
1398 ^^^^^
1400 .. mdp:: nwall
1402    (0)
1403    When set to 1 there is a wall at ``z=0``, when set to 2 there is
1404    also a wall at ``z=z-box``. Walls can only be used with :mdp:`pbc`
1405    ``=xy``. When set to 2, pressure coupling and Ewald summation can be
1406    used (it is usually best to use semiisotropic pressure coupling
1407    with the ``x/y`` compressibility set to 0, as otherwise the surface
1408    area will change). Walls interact wit the rest of the system
1409    through an optional :mdp:`wall-atomtype`. Energy groups ``wall0``
1410    and ``wall1`` (for :mdp:`nwall` =2) are added automatically to
1411    monitor the interaction of energy groups with each wall. The center
1412    of mass motion removal will be turned off in the ``z``-direction.
1414 .. mdp:: wall-atomtype
1416    the atom type name in the force field for each wall. By (for
1417    example) defining a special wall atom type in the topology with its
1418    own combination rules, this allows for independent tuning of the
1419    interaction of each atomtype with the walls.
1421 .. mdp:: wall-type
1423    .. mdp-value:: 9-3
1425       LJ integrated over the volume behind the wall: 9-3 potential
1427    .. mdp-value:: 10-4
1429       LJ integrated over the wall surface: 10-4 potential
1431    .. mdp-value:: 12-6
1433       direct LJ potential with the ``z`` distance from the wall
1435 .. mdp:: table
1437    user defined potentials indexed with the ``z`` distance from the
1438    wall, the tables are read analogously to the
1439    :mdp:`energygrp-table` option, where the first name is for a
1440    "normal" energy group and the second name is ``wall0`` or
1441    ``wall1``, only the dispersion and repulsion columns are used
1443 .. mdp:: wall-r-linpot
1445    (-1) [nm]
1446    Below this distance from the wall the potential is continued
1447    linearly and thus the force is constant. Setting this option to a
1448    postive value is especially useful for equilibration when some
1449    atoms are beyond a wall. When the value is <=0 (<0 for
1450    :mdp:`wall-type` =table), a fatal error is generated when atoms
1451    are beyond a wall.
1453 .. mdp:: wall-density
1455    [nm\ :sup:`-3`] / [nm\ :sup:`-2`]
1456    the number density of the atoms for each wall for wall types 9-3
1457    and 10-4
1459 .. mdp:: wall-ewald-zfac
1461    (3)
1462    The scaling factor for the third box vector for Ewald summation
1463    only, the minimum is 2. Ewald summation can only be used with
1464    :mdp:`nwall` =2, where one should use :mdp:`ewald-geometry`
1465    ``=3dc``. The empty layer in the box serves to decrease the
1466    unphysical Coulomb interaction between periodic images.
1469 COM pulling
1470 ^^^^^^^^^^^
1472 Note that where pulling coordinates are applicable, there can be more
1473 than one (set with :mdp:`pull-ncoords`) and multiple related :ref:`mdp`
1474 variables will exist accordingly. Documentation references to things
1475 like :mdp:`pull-coord1-vec` should be understood to apply to to the
1476 applicable pulling coordinate, eg. the second pull coordinate is described by
1477 pull-coord2-vec, pull-coord2-k, and so on.
1479 .. mdp:: pull
1481    .. mdp-value:: no
1483       No center of mass pulling. All the following pull options will
1484       be ignored (and if present in the :ref:`mdp` file, they unfortunately
1485       generate warnings)
1487    .. mdp-value:: yes
1489        Center of mass pulling will be applied on 1 or more groups using
1490        1 or more pull coordinates.
1492 .. mdp:: pull-cylinder-r
1494    (1.5) [nm]
1495    the radius of the cylinder for :mdp-value:`pull-coord1-geometry=cylinder`
1497 .. mdp:: pull-constr-tol
1499    (10\ :sup:`-6`)
1500    the relative constraint tolerance for constraint pulling
1502 .. mdp:: pull-print-com
1504    .. mdp-value:: no
1506       do not print the COM for any group
1508    .. mdp-value:: yes
1510       print the COM of all groups for all pull coordinates
1512 .. mdp:: pull-print-ref-value
1514    .. mdp-value:: no
1516       do not print the reference value for each pull coordinate
1518    .. mdp-value:: yes
1520       print the reference value for each pull coordinate
1522 .. mdp:: pull-print-components
1524    .. mdp-value:: no
1526       only print the distance for each pull coordinate
1528    .. mdp-value:: yes
1530       print the distance and Cartesian components selected in
1531       :mdp:`pull-coord1-dim`
1533 .. mdp:: pull-nstxout
1535    (50)
1536    frequency for writing out the COMs of all the pull group (0 is
1537    never)
1539 .. mdp:: pull-nstfout
1541    (50)
1542    frequency for writing out the force of all the pulled group
1543    (0 is never)
1545 .. mdp:: pull-pbc-ref-prev-step-com
1547    .. mdp-value:: no
1549       Use the reference atom (:mdp:`pull-group1-pbcatom`) for the
1550       treatment of periodic boundary conditions.
1552    .. mdp-value:: yes
1554       Use the COM of the previous step as reference for the treatment
1555       of periodic boundary conditions. The reference is initialized
1556       using the reference atom (:mdp:`pull-group1-pbcatom`), which should
1557       be located centrally in the group. Using the COM from the
1558       previous step can be useful if one or more pull groups are large.
1560 .. mdp:: pull-xout-average
1562    .. mdp-value:: no
1564       Write the instantaneous coordinates for all the pulled groups.
1566    .. mdp-value:: yes
1568       Write the average coordinates (since last output) for all the
1569       pulled groups. N.b., some analysis tools might expect instantaneous
1570       pull output.
1572 .. mdp:: pull-fout-average
1574    .. mdp-value:: no
1576       Write the instantaneous force for all the pulled groups.
1578    .. mdp-value:: yes
1580       Write the average force (since last output) for all the
1581       pulled groups. N.b., some analysis tools might expect instantaneous
1582       pull output.
1584 .. mdp:: pull-ngroups
1586    (1)
1587    The number of pull groups, not including the absolute reference
1588    group, when used. Pull groups can be reused in multiple pull
1589    coordinates. Below only the pull options for group 1 are given,
1590    further groups simply increase the group index number.
1592 .. mdp:: pull-ncoords
1594    (1)
1595    The number of pull coordinates. Below only the pull options for
1596    coordinate 1 are given, further coordinates simply increase the
1597    coordinate index number.
1599 .. mdp:: pull-group1-name
1601    The name of the pull group, is looked up in the index file or in
1602    the default groups to obtain the atoms involved.
1604 .. mdp:: pull-group1-weights
1606    Optional relative weights which are multiplied with the masses of
1607    the atoms to give the total weight for the COM. The number should
1608    be 0, meaning all 1, or the number of atoms in the pull group.
1610 .. mdp:: pull-group1-pbcatom
1612    (0)
1613    The reference atom for the treatment of periodic boundary
1614    conditions inside the group (this has no effect on the treatment of
1615    the pbc between groups). This option is only important when the
1616    diameter of the pull group is larger than half the shortest box
1617    vector. For determining the COM, all atoms in the group are put at
1618    their periodic image which is closest to
1619    :mdp:`pull-group1-pbcatom`. A value of 0 means that the middle
1620    atom (number wise) is used, which is only safe for small groups.
1621    :ref:`gmx grompp` checks that the maximum distance from the reference
1622    atom (specifically chosen, or not) to the other atoms in the group
1623    is not too large. This parameter is not used with
1624    :mdp:`pull-coord1-geometry` cylinder. A value of -1 turns on cosine
1625    weighting, which is useful for a group of molecules in a periodic
1626    system, *e.g.* a water slab (see Engin et al. J. Chem. Phys. B
1627    2010).
1629 .. mdp:: pull-coord1-type
1631    .. mdp-value:: umbrella
1633       Center of mass pulling using an umbrella potential between the
1634       reference group and one or more groups.
1636    .. mdp-value:: constraint
1638       Center of mass pulling using a constraint between the reference
1639       group and one or more groups. The setup is identical to the
1640       option umbrella, except for the fact that a rigid constraint is
1641       applied instead of a harmonic potential.
1643    .. mdp-value:: constant-force
1645       Center of mass pulling using a linear potential and therefore a
1646       constant force. For this option there is no reference position
1647       and therefore the parameters :mdp:`pull-coord1-init` and
1648       :mdp:`pull-coord1-rate` are not used.
1650    .. mdp-value:: flat-bottom
1652       At distances above :mdp:`pull-coord1-init` a harmonic potential
1653       is applied, otherwise no potential is applied.
1655    .. mdp-value:: flat-bottom-high
1657       At distances below :mdp:`pull-coord1-init` a harmonic potential
1658       is applied, otherwise no potential is applied.
1660    .. mdp-value:: external-potential
1662       An external potential that needs to be provided by another
1663       module.
1665 .. mdp:: pull-coord1-potential-provider
1667       The name of the external module that provides the potential for
1668       the case where :mdp:`pull-coord1-type` is external-potential.
1670 .. mdp:: pull-coord1-geometry
1672    .. mdp-value:: distance
1674       Pull along the vector connecting the two groups. Components can
1675       be selected with :mdp:`pull-coord1-dim`.
1677    .. mdp-value:: direction
1679       Pull in the direction of :mdp:`pull-coord1-vec`.
1681    .. mdp-value:: direction-periodic
1683       As :mdp-value:`pull-coord1-geometry=direction`, but does not apply
1684       periodic box vector corrections to keep the distance within half
1685       the box length. This is (only) useful for pushing groups apart
1686       by more than half the box length by continuously changing the reference
1687       location using a pull rate. With this geometry the box should not be
1688       dynamic (*e.g.* no pressure scaling) in the pull dimensions and
1689       the pull force is not added to the virial.
1691    .. mdp-value:: direction-relative
1693       As :mdp-value:`pull-coord1-geometry=direction`, but the pull vector is the vector
1694       that points from the COM of a third to the COM of a fourth pull
1695       group. This means that 4 groups need to be supplied in
1696       :mdp:`pull-coord1-groups`. Note that the pull force will give
1697       rise to a torque on the pull vector, which is turn leads to
1698       forces perpendicular to the pull vector on the two groups
1699       defining the vector. If you want a pull group to move between
1700       the two groups defining the vector, simply use the union of
1701       these two groups as the reference group.
1703    .. mdp-value:: cylinder
1705       Designed for pulling with respect to a layer where the reference
1706       COM is given by a local cylindrical part of the reference group.
1707       The pulling is in the direction of :mdp:`pull-coord1-vec`. From
1708       the first of the two groups in :mdp:`pull-coord1-groups` a
1709       cylinder is selected around the axis going through the COM of
1710       the second group with direction :mdp:`pull-coord1-vec` with
1711       radius :mdp:`pull-cylinder-r`. Weights of the atoms decrease
1712       continously to zero as the radial distance goes from 0 to
1713       :mdp:`pull-cylinder-r` (mass weighting is also used). The radial
1714       dependence gives rise to radial forces on both pull groups.
1715       Note that the radius should be smaller than half the box size.
1716       For tilted cylinders they should be even smaller than half the
1717       box size since the distance of an atom in the reference group
1718       from the COM of the pull group has both a radial and an axial
1719       component. This geometry is not supported with constraint
1720       pulling.
1722    .. mdp-value:: angle
1724       Pull along an angle defined by four groups. The angle is
1725       defined as the angle between two vectors: the vector connecting
1726       the COM of the first group to the COM of the second group and
1727       the vector connecting the COM of the third group to the COM of
1728       the fourth group.
1730    .. mdp-value:: angle-axis
1732       As :mdp-value:`pull-coord1-geometry=angle` but the second vector is given by :mdp:`pull-coord1-vec`.
1733       Thus, only the two groups that define the first vector need to be given.
1735    .. mdp-value:: dihedral
1737       Pull along a dihedral angle defined by six groups. These pairwise
1738       define three vectors: the vector connecting the COM of group 1
1739       to the COM of group 2, the COM of group 3 to the COM of group 4,
1740       and the COM of group 5 to the COM group 6. The dihedral angle is
1741       then defined as the angle between two planes: the plane spanned by the
1742       the two first vectors and the plane spanned the two last vectors.
1745 .. mdp:: pull-coord1-groups
1747    The group indices on which this pull coordinate will operate.
1748    The number of group indices required is geometry dependent.
1749    The first index can be 0, in which case an
1750    absolute reference of :mdp:`pull-coord1-origin` is used. With an
1751    absolute reference the system is no longer translation invariant
1752    and one should think about what to do with the center of mass
1753    motion.
1755 .. mdp:: pull-coord1-dim
1757    (Y Y Y)
1758    Selects the dimensions that this pull coordinate acts on and that
1759    are printed to the output files when
1760    :mdp:`pull-print-components` = :mdp-value:`pull-coord1-start=yes`. With
1761    :mdp:`pull-coord1-geometry` = :mdp-value:`pull-coord1-geometry=distance`, only Cartesian
1762    components set to Y contribute to the distance. Thus setting this
1763    to Y Y N results in a distance in the x/y plane. With other
1764    geometries all dimensions with non-zero entries in
1765    :mdp:`pull-coord1-vec` should be set to Y, the values for other
1766    dimensions only affect the output.
1768 .. mdp:: pull-coord1-origin
1770    (0.0 0.0 0.0)
1771    The pull reference position for use with an absolute reference.
1773 .. mdp:: pull-coord1-vec
1775    (0.0 0.0 0.0)
1776    The pull direction. :ref:`gmx grompp` normalizes the vector.
1778 .. mdp:: pull-coord1-start
1780    .. mdp-value:: no
1782       do not modify :mdp:`pull-coord1-init`
1784    .. mdp-value:: yes
1786       add the COM distance of the starting conformation to
1787       :mdp:`pull-coord1-init`
1789 .. mdp:: pull-coord1-init
1791    (0.0) [nm] or [deg]
1792    The reference distance or reference angle at t=0.
1794 .. mdp:: pull-coord1-rate
1796    (0) [nm/ps] or [deg/ps]
1797    The rate of change of the reference position or reference angle.
1799 .. mdp:: pull-coord1-k
1801    (0) [kJ mol\ :sup:`-1` nm\ :sup:`-2`] or [kJ mol\ :sup:`-1` nm\ :sup:`-1`] or
1802    [kJ mol\ :sup:`-1` rad\ :sup:`-2`] or [kJ mol\ :sup:`-1` rad\ :sup:`-1`]
1803    The force constant. For umbrella pulling this is the harmonic force
1804    constant in kJ mol\ :sup:`-1` nm\ :sup:`-2` (or kJ mol\ :sup:`-1` rad\ :sup:`-2`
1805    for angles). For constant force pulling this is the
1806    force constant of the linear potential, and thus the negative (!)
1807    of the constant force in kJ mol\ :sup:`-1` nm\ :sup:`-1`
1808    (or kJ mol\ :sup:`-1` rad\ :sup:`-1` for angles).
1809    Note that for angles the force constant is expressed in terms of radians
1810    (while :mdp:`pull-coord1-init` and :mdp:`pull-coord1-rate` are expressed in degrees).
1812 .. mdp:: pull-coord1-kB
1814    (pull-k1) [kJ mol\ :sup:`-1` nm\ :sup:`-2`] or [kJ mol\ :sup:`-1` nm\ :sup:`-1`]
1815    or [kJ mol\ :sup:`-1` rad\ :sup:`-2`] or [kJ mol\ :sup:`-1` rad\ :sup:`-1`]
1816    As :mdp:`pull-coord1-k`, but for state B. This is only used when
1817    :mdp:`free-energy` is turned on. The force constant is then (1 -
1818    lambda) * :mdp:`pull-coord1-k` + lambda * :mdp:`pull-coord1-kB`.
1820 AWH adaptive biasing
1821 ^^^^^^^^^^^^^^^^^^^^
1823 .. mdp:: awh
1825    .. mdp-value:: no
1827       No biasing.
1829    .. mdp-value:: yes
1831       Adaptively bias a reaction coordinate using the AWH method and estimate
1832       the corresponding PMF. The PMF and other AWH data are written to energy
1833       file at an interval set by :mdp:`awh-nstout` and can be extracted with
1834       the ``gmx awh`` tool. The AWH coordinate can be
1835       multidimensional and is defined by mapping each dimension to a pull coordinate index.
1836       This is only allowed if :mdp-value:`pull-coord1-type=external-potential` and
1837       :mdp:`pull-coord1-potential-provider` = ``awh`` for the concerned pull coordinate
1838       indices. Pull geometry 'direction-periodic' is not supported by AWH.
1840 .. mdp:: awh-potential
1842    .. mdp-value:: convolved
1844       The applied biasing potential is the convolution of the bias function and a
1845       set of harmonic umbrella potentials (see :mdp-value:`awh-potential=umbrella` below). This results
1846       in a smooth potential function and force. The resolution of the potential is set
1847       by the force constant of each umbrella, see :mdp:`awh1-dim1-force-constant`.
1849    .. mdp-value:: umbrella
1851       The potential bias is applied by controlling the position of an harmonic potential
1852       using Monte-Carlo sampling.  The force constant is set with
1853       :mdp:`awh1-dim1-force-constant`. The umbrella location
1854       is sampled using Monte-Carlo every :mdp:`awh-nstsample` steps.
1855       There are no advantages to using an umbrella.
1856       This option is mainly for comparison and testing purposes.
1858 .. mdp:: awh-share-multisim
1860    .. mdp-value:: no
1862       AWH will not share biases across simulations started with
1863       :ref:`gmx mdrun` option ``-multidir``. The biases will be independent.
1865    .. mdp-value:: yes
1867       With :ref:`gmx mdrun` and option ``-multidir`` the bias and PMF estimates
1868       for biases with :mdp:`awh1-share-group` >0 will be shared across simulations
1869       with the biases with the same :mdp:`awh1-share-group` value.
1870       The simulations should have the same AWH settings for sharing to make sense.
1871       :ref:`gmx mdrun` will check whether the simulations are technically
1872       compatible for sharing, but the user should check that bias sharing
1873       physically makes sense.
1875 .. mdp:: awh-seed
1877    (-1) Random seed for Monte-Carlo sampling the umbrella position,
1878    where -1 indicates to generate a seed. Only used with
1879    :mdp-value:`awh-potential=umbrella`.
1881 .. mdp:: awh-nstout
1883    (100000)
1884    Number of steps between printing AWH data to the energy file, should be
1885    a multiple of :mdp:`nstenergy`.
1887 .. mdp:: awh-nstsample
1889    (10)
1890    Number of steps between sampling of the coordinate value. This sampling
1891    is the basis for updating the bias and estimating the PMF and other AWH observables.
1893 .. mdp:: awh-nsamples-update
1895    (10)
1896    The number of coordinate samples used for each AWH update.
1897    The update interval in steps is :mdp:`awh-nstsample` times this value.
1899 .. mdp:: awh-nbias
1901    (1)
1902    The number of biases, each acting on its own coordinate.
1903    The following options should be specified
1904    for each bias although below only the options for bias number 1 is shown. Options for
1905    other bias indices are  obtained by replacing '1' by the bias index.
1907 .. mdp:: awh1-error-init
1909    (10.0) [kJ mol\ :sup:`-1`]
1910    Estimated initial average error of the PMF for this bias. This value together with the
1911    given diffusion constant(s) :mdp:`awh1-dim1-diffusion` determine the initial biasing rate.
1912    The error is obviously not known *a priori*. Only a rough estimate of :mdp:`awh1-error-init`
1913    is needed however.
1914    As a  general guideline, leave :mdp:`awh1-error-init` to its default value when starting a new
1915    simulation. On the other hand, when there is *a priori* knowledge of the PMF (e.g. when
1916    an initial PMF estimate is provided, see the :mdp:`awh1-user-data` option)
1917    then :mdp:`awh1-error-init` should reflect that knowledge.
1919 .. mdp:: awh1-growth
1921    .. mdp-value:: exp-linear
1923    Each bias keeps a reference weight histogram for the coordinate samples.
1924    Its size sets the magnitude of the bias function and free energy estimate updates
1925    (few samples corresponds to large updates and vice versa).
1926    Thus, its growth rate sets the maximum convergence rate.
1927    By default, there is an initial stage in which the histogram grows close to exponentially (but slower than the sampling rate).
1928    In the final stage that follows, the growth rate is linear and equal to the sampling rate (set by :mdp:`awh-nstsample`).
1929    The initial stage is typically necessary for efficient convergence when starting a new simulation where
1930    high free energy barriers have not yet been flattened by the bias.
1932    .. mdp-value:: linear
1934    As :mdp-value:`awh1-growth=exp-linear` but skip the initial stage. This may be useful if there is *a priori*
1935    knowledge (see :mdp:`awh1-error-init`) which eliminates the need for an initial stage. This is also
1936    the setting compatible with :mdp-value:`awh1-target=local-boltzmann`.
1938 .. mdp:: awh1-equilibrate-histogram
1940    .. mdp-value:: no
1942       Do not equilibrate histogram.
1944    .. mdp-value:: yes
1946       Before entering the initial stage (see :mdp-value:`awh1-growth=exp-linear`), make sure the
1947       histogram of sampled weights is following the target distribution closely enough (specifically,
1948       at least 80% of the target region needs to have a local relative error of less than 20%). This
1949       option would typically only be used when :mdp:`awh1-share-group` > 0
1950       and the initial configurations poorly represent the target
1951       distribution.
1953 .. mdp:: awh1-target
1955    .. mdp-value:: constant
1957       The bias is tuned towards a constant (uniform) coordinate distribution
1958       in the defined sampling interval (defined by  [:mdp:`awh1-dim1-start`, :mdp:`awh1-dim1-end`]).
1960    .. mdp-value:: cutoff
1962       Similar to :mdp-value:`awh1-target=constant`, but the target
1963       distribution is proportional to 1/(1 + exp(F - :mdp-value:`awh1-target=cutoff`)),
1964       where F is the free energy relative to the estimated global minimum.
1965       This provides a smooth switch of a flat target distribution in
1966       regions with free energy lower than the cut-off to a Boltzmann
1967       distribution in regions with free energy higher than the cut-off.
1969    .. mdp-value:: boltzmann
1971       The target distribution is a Boltzmann distribtution with a scaled beta (inverse temperature)
1972       factor given by :mdp:`awh1-target-beta-scaling`. *E.g.*, a value of 0.1
1973       would give the same coordinate distribution as sampling with a simulation temperature
1974       scaled by 10.
1976    .. mdp-value:: local-boltzmann
1978       Same target distribution and use of :mdp:`awh1-target-beta-scaling`
1979       but the convergence towards the target distribution is inherently local *i.e.*, the rate of
1980       change of the bias only depends on the local sampling. This local convergence property is
1981       only compatible with :mdp-value:`awh1-growth=linear`, since for
1982       :mdp-value:`awh1-growth=exp-linear` histograms are globally rescaled in the initial stage.
1984 .. mdp:: awh1-target-beta-scaling
1986    (0)
1987    For :mdp-value:`awh1-target=boltzmann` and :mdp-value:`awh1-target=local-boltzmann`
1988    it is the unitless beta scaling factor taking values in (0,1).
1990 .. mdp:: awh1-target-cutoff
1992    (0) [kJ mol\ :sup:`-1`]
1993    For :mdp-value:`awh1-target=cutoff` this is the cutoff, should be > 0.
1995 .. mdp:: awh1-user-data
1997    .. mdp-value:: no
1999       Initialize the PMF and target distribution with default values.
2001    .. mdp-value:: yes
2003       Initialize the PMF and target distribution with user provided data. For :mdp:`awh-nbias` = 1,
2004       :ref:`gmx mdrun` will expect a file ``awhinit.xvg`` to be present in the run directory.
2005       For multiple biases, :ref:`gmx mdrun` expects files ``awhinit1.xvg``, ``awhinit2.xvg``, etc.
2006       The file name can be changed with the ``-awh`` option.
2007       The first :mdp:`awh1-ndim` columns of
2008       each input file should contain the coordinate values, such that each row defines a point in
2009       coordinate space. Column :mdp:`awh1-ndim` + 1 should contain the PMF value (in kT) for each point.
2010       The target distribution column can either follow the PMF (column  :mdp:`awh1-ndim` + 2) or
2011       be in the same column as written by :ref:`gmx awh`.
2013 .. mdp:: awh1-share-group
2015    .. mdp-value:: 0
2017       Do not share the bias.
2019    .. mdp-value:: positive
2021       Share the bias and PMF estimates within and/or between simulations.
2022       Within a simulation, the bias will be shared between biases that have the
2023       same :mdp:`awh1-share-group` index (note that the current code does not support this).
2024       With :mdp-value:`awh-share-multisim=yes` and
2025       :ref:`gmx mdrun` option ``-multidir`` the bias will also be shared across simulations.
2026       Sharing may increase convergence initially, although the starting configurations
2027       can be critical, especially when sharing between many biases.
2028       Currently, positive group values should start at 1 and increase
2029       by 1 for each subsequent bias that is shared.
2031 .. mdp:: awh1-ndim
2033    (1) [integer]
2034    Number of dimensions of the coordinate, each dimension maps to 1 pull coordinate.
2035    The following options should be specified for each such dimension. Below only
2036    the options for dimension number 1 is shown. Options for other dimension indices are
2037    obtained by replacing '1' by the dimension index.
2039 .. mdp:: awh1-dim1-coord-provider
2041    .. mdp-value:: pull
2043       The module providing the reaction coordinate for this dimension.
2044       Currently AWH can only act on pull coordinates.
2046 .. mdp:: awh1-dim1-coord-index
2048    (1)
2049    Index of the pull coordinate defining this coordinate dimension.
2051 .. mdp:: awh1-dim1-force-constant
2053    (0) [kJ mol\ :sup:`-1` nm\ :sup:`-2`] or [kJ mol\ :sup:`-1` rad\ :sup:`-2`]
2054    Force constant for the (convolved) umbrella potential(s) along this
2055    coordinate dimension.
2057 .. mdp:: awh1-dim1-start
2059    (0.0) [nm] or [rad]
2060    Start value of the sampling interval along this dimension. The range of allowed
2061    values depends on the relevant pull geometry (see :mdp:`pull-coord1-geometry`).
2062    For dihedral geometries :mdp:`awh1-dim1-start` greater than :mdp:`awh1-dim1-end`
2063    is allowed. The interval will then wrap around from +period/2 to -period/2.
2064    For the direction geometry, the dimension is made periodic when
2065    the direction is along a box vector and covers more than 95%
2066    of the box length. Note that one should not apply pressure coupling
2067    along a periodic dimension.
2069 .. mdp:: awh1-dim1-end
2071    (0.0) [nm] or [rad]
2072    End value defining the sampling interval together with :mdp:`awh1-dim1-start`.
2074 .. mdp:: awh1-dim1-diffusion
2076    (10\ :sup:`-5`) [nm\ :sup:`2`/ps] or [rad\ :sup:`2`/ps]
2077    Estimated diffusion constant for this coordinate dimension determining the initial
2078    biasing rate. This needs only be a rough estimate and should not critically
2079    affect the results unless it is set to something very low, leading to slow convergence,
2080    or very high, forcing the system far from equilibrium. Not setting this value
2081    explicitly generates a warning.
2083 .. mdp:: awh1-dim1-cover-diameter
2085    (0.0) [nm] or [rad]
2086    Diameter that needs to be sampled by a single simulation around a coordinate value
2087    before the point is considered covered in the initial stage (see :mdp-value:`awh1-growth=exp-linear`).
2088    A value > 0  ensures that for each covering there is a continuous transition of this diameter
2089    across each coordinate value.
2090    This is trivially true for independent simulations but not for for multiple bias-sharing simulations
2091    (:mdp:`awh1-share-group`>0).
2092    For a diameter = 0, covering occurs as soon as the simulations have sampled the whole interval, which
2093    for many sharing simulations does not guarantee transitions across free energy barriers.
2094    On the other hand, when the diameter >= the sampling interval length, covering occurs when a single simulation
2095    has independently sampled the whole interval.
2097 Enforced rotation
2098 ^^^^^^^^^^^^^^^^^
2100 These :ref:`mdp` parameters can be used enforce the rotation of a group of atoms,
2101 e.g. a protein subunit. The `reference manual`_ describes in detail 13 different potentials
2102 that can be used to achieve such a rotation.
2104 .. mdp:: rotation
2106    .. mdp-value:: no
2108       No enforced rotation will be applied. All enforced rotation options will
2109       be ignored (and if present in the :ref:`mdp` file, they unfortunately
2110       generate warnings).
2112    .. mdp-value:: yes
2114       Apply the rotation potential specified by :mdp:`rot-type0` to the group of atoms given
2115       under the :mdp:`rot-group0` option.
2117 .. mdp:: rot-ngroups
2119    (1)
2120    Number of rotation groups.
2122 .. mdp:: rot-group0
2124    Name of rotation group 0 in the index file.
2126 .. mdp:: rot-type0
2128    (iso)
2129    Type of rotation potential that is applied to rotation group 0. Can be of of the following:
2130    ``iso``, ``iso-pf``, ``pm``, ``pm-pf``, ``rm``, ``rm-pf``, ``rm2``, ``rm2-pf``,
2131    ``flex``, ``flex-t``, ``flex2``, or ``flex2-t``.
2133 .. mdp:: rot-massw0
2135    (no)
2136    Use mass weighted rotation group positions.
2138 .. mdp:: rot-vec0
2140    (1.0 0.0 0.0)
2141    Rotation vector, will get normalized.
2143 .. mdp:: rot-pivot0
2145    (0.0 0.0 0.0) [nm]
2146    Pivot point for the potentials ``iso``, ``pm``, ``rm``, and ``rm2``.
2148 .. mdp:: rot-rate0
2150    (0) [degree ps\ :sup:`-1`]
2151    Reference rotation rate of group 0.
2153 .. mdp:: rot-k0
2155    (0) [kJ mol\ :sup:`-1` nm\ :sup:`-2`]
2156    Force constant for group 0.
2158 .. mdp:: rot-slab-dist0
2160    (1.5) [nm]
2161    Slab distance, if a flexible axis rotation type was chosen.
2163 .. mdp:: rot-min-gauss0
2165    (0.001)
2166    Minimum value (cutoff) of Gaussian function for the force to be evaluated
2167    (for the flexible axis potentials).
2169 .. mdp:: rot-eps0
2171    (0.0001) [nm\ :sup:`2`]
2172    Value of additive constant epsilon for ``rm2*`` and ``flex2*`` potentials.
2174 .. mdp:: rot-fit-method0
2176    (rmsd)
2177    Fitting method when determining the actual angle of a rotation group
2178    (can be one of ``rmsd``, ``norm``, or ``potential``).
2180 .. mdp:: rot-potfit-nsteps0
2182    (21)
2183    For fit type ``potential``, the number of angular positions around the reference angle for which the
2184    rotation potential is evaluated.
2186 .. mdp:: rot-potfit-step0
2188    (0.25)
2189    For fit type ``potential``, the distance in degrees between two angular positions.
2191 .. mdp:: rot-nstrout
2193    (100)
2194    Output frequency (in steps) for the angle of the rotation group, as well as for the torque
2195    and the rotation potential energy.
2197 .. mdp:: rot-nstsout
2199    (1000)
2200    Output frequency for per-slab data of the flexible axis potentials, i.e. angles, torques and slab centers.
2203 NMR refinement
2204 ^^^^^^^^^^^^^^
2206 .. mdp:: disre
2208    .. mdp-value:: no
2210       ignore distance restraint information in topology file
2212    .. mdp-value:: simple
2214       simple (per-molecule) distance restraints.
2216    .. mdp-value:: ensemble
2218       distance restraints over an ensemble of molecules in one
2219       simulation box. Normally, one would perform ensemble averaging
2220       over multiple simulations, using ``mdrun
2221       -multidir``. The environment
2222       variable ``GMX_DISRE_ENSEMBLE_SIZE`` sets the number of systems
2223       within each ensemble (usually equal to the number of directories
2224       supplied to ``mdrun -multidir``).
2226 .. mdp:: disre-weighting
2228    .. mdp-value:: equal
2230       divide the restraint force equally over all atom pairs in the
2231       restraint
2233    .. mdp-value:: conservative
2235       the forces are the derivative of the restraint potential, this
2236       results in an weighting of the atom pairs to the reciprocal
2237       seventh power of the displacement. The forces are conservative
2238       when :mdp:`disre-tau` is zero.
2240 .. mdp:: disre-mixed
2242    .. mdp-value:: no
2244       the violation used in the calculation of the restraint force is
2245       the time-averaged violation
2247    .. mdp-value:: yes
2249       the violation used in the calculation of the restraint force is
2250       the square root of the product of the time-averaged violation
2251       and the instantaneous violation
2253 .. mdp:: disre-fc
2255    (1000) [kJ mol\ :sup:`-1` nm\ :sup:`-2`]
2256    force constant for distance restraints, which is multiplied by a
2257    (possibly) different factor for each restraint given in the ``fac``
2258    column of the interaction in the topology file.
2260 .. mdp:: disre-tau
2262    (0) [ps]
2263    time constant for distance restraints running average. A value of
2264    zero turns off time averaging.
2266 .. mdp:: nstdisreout
2268    (100) [steps]
2269    period between steps when the running time-averaged and
2270    instantaneous distances of all atom pairs involved in restraints
2271    are written to the energy file (can make the energy file very
2272    large)
2274 .. mdp:: orire
2276    .. mdp-value:: no
2278       ignore orientation restraint information in topology file
2280    .. mdp-value:: yes
2282       use orientation restraints, ensemble averaging can be performed
2283       with ``mdrun -multidir``
2285 .. mdp:: orire-fc
2287    (0) [kJ mol\ :sup:`-1`]
2288    force constant for orientation restraints, which is multiplied by a
2289    (possibly) different weight factor for each restraint, can be set
2290    to zero to obtain the orientations from a free simulation
2292 .. mdp:: orire-tau
2294    (0) [ps]
2295    time constant for orientation restraints running average. A value
2296    of zero turns off time averaging.
2298 .. mdp:: orire-fitgrp
2300    fit group for orientation restraining. This group of atoms is used
2301    to determine the rotation **R** of the system with respect to the
2302    reference orientation. The reference orientation is the starting
2303    conformation of the first subsystem. For a protein, backbone is a
2304    reasonable choice
2306 .. mdp:: nstorireout
2308    (100) [steps]
2309    period between steps when the running time-averaged and
2310    instantaneous orientations for all restraints, and the molecular
2311    order tensor are written to the energy file (can make the energy
2312    file very large)
2315 Free energy calculations
2316 ^^^^^^^^^^^^^^^^^^^^^^^^
2318 .. mdp:: free-energy
2320    .. mdp-value:: no
2322       Only use topology A.
2324    .. mdp-value:: yes
2326       Interpolate between topology A (lambda=0) to topology B
2327       (lambda=1) and write the derivative of the Hamiltonian with
2328       respect to lambda (as specified with :mdp:`dhdl-derivatives`),
2329       or the Hamiltonian differences with respect to other lambda
2330       values (as specified with foreign lambda) to the energy file
2331       and/or to ``dhdl.xvg``, where they can be processed by, for
2332       example :ref:`gmx bar`. The potentials, bond-lengths and angles
2333       are interpolated linearly as described in the manual. When
2334       :mdp:`sc-alpha` is larger than zero, soft-core potentials are
2335       used for the LJ and Coulomb interactions.
2337 .. mdp:: expanded
2339    Turns on expanded ensemble simulation, where the alchemical state
2340    becomes a dynamic variable, allowing jumping between different
2341    Hamiltonians. See the expanded ensemble options for controlling how
2342    expanded ensemble simulations are performed. The different
2343    Hamiltonians used in expanded ensemble simulations are defined by
2344    the other free energy options.
2346 .. mdp:: init-lambda
2348    (-1)
2349    starting value for lambda (float). Generally, this should only be
2350    used with slow growth (*i.e.* nonzero :mdp:`delta-lambda`). In
2351    other cases, :mdp:`init-lambda-state` should be specified
2352    instead. Must be greater than or equal to 0.
2354 .. mdp:: delta-lambda
2356    (0)
2357    increment per time step for lambda
2359 .. mdp:: init-lambda-state
2361    (-1)
2362    starting value for the lambda state (integer). Specifies which
2363    columm of the lambda vector (:mdp:`coul-lambdas`,
2364    :mdp:`vdw-lambdas`, :mdp:`bonded-lambdas`,
2365    :mdp:`restraint-lambdas`, :mdp:`mass-lambdas`,
2366    :mdp:`temperature-lambdas`, :mdp:`fep-lambdas`) should be
2367    used. This is a zero-based index: :mdp:`init-lambda-state` 0 means
2368    the first column, and so on.
2370 .. mdp:: fep-lambdas
2372    [array]
2373    Zero, one or more lambda values for which Delta H values will be
2374    determined and written to dhdl.xvg every :mdp:`nstdhdl`
2375    steps. Values must be between 0 and 1. Free energy differences
2376    between different lambda values can then be determined with
2377    :ref:`gmx bar`. :mdp:`fep-lambdas` is different from the
2378    other -lambdas keywords because all components of the lambda vector
2379    that are not specified will use :mdp:`fep-lambdas` (including
2380    :mdp:`restraint-lambdas` and therefore the pull code restraints).
2382 .. mdp:: coul-lambdas
2384    [array]
2385    Zero, one or more lambda values for which Delta H values will be
2386    determined and written to dhdl.xvg every :mdp:`nstdhdl`
2387    steps. Values must be between 0 and 1. Only the electrostatic
2388    interactions are controlled with this component of the lambda
2389    vector (and only if the lambda=0 and lambda=1 states have differing
2390    electrostatic interactions).
2392 .. mdp:: vdw-lambdas
2394    [array]
2395    Zero, one or more lambda values for which Delta H values will be
2396    determined and written to dhdl.xvg every :mdp:`nstdhdl`
2397    steps. Values must be between 0 and 1. Only the van der Waals
2398    interactions are controlled with this component of the lambda
2399    vector.
2401 .. mdp:: bonded-lambdas
2403    [array]
2404    Zero, one or more lambda values for which Delta H values will be
2405    determined and written to dhdl.xvg every :mdp:`nstdhdl`
2406    steps. Values must be between 0 and 1. Only the bonded interactions
2407    are controlled with this component of the lambda vector.
2409 .. mdp:: restraint-lambdas
2411    [array]
2412    Zero, one or more lambda values for which Delta H values will be
2413    determined and written to dhdl.xvg every :mdp:`nstdhdl`
2414    steps. Values must be between 0 and 1. Only the restraint
2415    interactions: dihedral restraints, and the pull code restraints are
2416    controlled with this component of the lambda vector.
2418 .. mdp:: mass-lambdas
2420    [array]
2421    Zero, one or more lambda values for which Delta H values will be
2422    determined and written to dhdl.xvg every :mdp:`nstdhdl`
2423    steps. Values must be between 0 and 1. Only the particle masses are
2424    controlled with this component of the lambda vector.
2426 .. mdp:: temperature-lambdas
2428    [array]
2429    Zero, one or more lambda values for which Delta H values will be
2430    determined and written to dhdl.xvg every :mdp:`nstdhdl`
2431    steps. Values must be between 0 and 1. Only the temperatures
2432    controlled with this component of the lambda vector. Note that
2433    these lambdas should not be used for replica exchange, only for
2434    simulated tempering.
2436 .. mdp:: calc-lambda-neighbors
2438    (1)
2439    Controls the number of lambda values for which Delta H values will
2440    be calculated and written out, if :mdp:`init-lambda-state` has
2441    been set. A positive value will limit the number of lambda points
2442    calculated to only the nth neighbors of :mdp:`init-lambda-state`:
2443    for example, if :mdp:`init-lambda-state` is 5 and this parameter
2444    has a value of 2, energies for lambda points 3-7 will be calculated
2445    and writen out. A value of -1 means all lambda points will be
2446    written out. For normal BAR such as with :ref:`gmx bar`, a value of
2447    1 is sufficient, while for MBAR -1 should be used.
2449 .. mdp:: sc-alpha
2451    (0)
2452    the soft-core alpha parameter, a value of 0 results in linear
2453    interpolation of the LJ and Coulomb interactions
2455 .. mdp:: sc-r-power
2457    (6)
2458    power 6 for the radial term in the soft-core equation.
2460 .. mdp:: sc-coul
2462    (no)
2463    Whether to apply the soft-core free energy interaction
2464    transformation to the Columbic interaction of a molecule. Default
2465    is no, as it is generally more efficient to turn off the Coulomic
2466    interactions linearly before turning off the van der Waals
2467    interactions. Note that it is only taken into account when lambda
2468    states are used, not with :mdp:`couple-lambda0` /
2469    :mdp:`couple-lambda1`, and you can still turn off soft-core
2470    interactions by setting :mdp:`sc-alpha` to 0.
2472 .. mdp:: sc-power
2474    (0)
2475    the power for lambda in the soft-core function, only the values 1
2476    and 2 are supported
2478 .. mdp:: sc-sigma
2480    (0.3) [nm]
2481    the soft-core sigma for particles which have a C6 or C12 parameter
2482    equal to zero or a sigma smaller than :mdp:`sc-sigma`
2484 .. mdp:: couple-moltype
2486    Here one can supply a molecule type (as defined in the topology)
2487    for calculating solvation or coupling free energies. There is a
2488    special option ``system`` that couples all molecule types in the
2489    system. This can be useful for equilibrating a system starting from
2490    (nearly) random coordinates. :mdp:`free-energy` has to be turned
2491    on. The Van der Waals interactions and/or charges in this molecule
2492    type can be turned on or off between lambda=0 and lambda=1,
2493    depending on the settings of :mdp:`couple-lambda0` and
2494    :mdp:`couple-lambda1`. If you want to decouple one of several
2495    copies of a molecule, you need to copy and rename the molecule
2496    definition in the topology.
2498 .. mdp:: couple-lambda0
2500    .. mdp-value:: vdw-q
2502       all interactions are on at lambda=0
2504    .. mdp-value:: vdw
2506       the charges are zero (no Coulomb interactions) at lambda=0
2508    .. mdp-value:: q
2510       the Van der Waals interactions are turned at lambda=0; soft-core
2511       interactions will be required to avoid singularities
2513    .. mdp-value:: none
2515       the Van der Waals interactions are turned off and the charges
2516       are zero at lambda=0; soft-core interactions will be required to
2517       avoid singularities.
2519 .. mdp:: couple-lambda1
2521    analogous to :mdp:`couple-lambda1`, but for lambda=1
2523 .. mdp:: couple-intramol
2525    .. mdp-value:: no
2527       All intra-molecular non-bonded interactions for moleculetype
2528       :mdp:`couple-moltype` are replaced by exclusions and explicit
2529       pair interactions. In this manner the decoupled state of the
2530       molecule corresponds to the proper vacuum state without
2531       periodicity effects.
2533    .. mdp-value:: yes
2535       The intra-molecular Van der Waals and Coulomb interactions are
2536       also turned on/off. This can be useful for partitioning
2537       free-energies of relatively large molecules, where the
2538       intra-molecular non-bonded interactions might lead to
2539       kinetically trapped vacuum conformations. The 1-4 pair
2540       interactions are not turned off.
2542 .. mdp:: nstdhdl
2544    (100)
2545    the frequency for writing dH/dlambda and possibly Delta H to
2546    dhdl.xvg, 0 means no ouput, should be a multiple of
2547    :mdp:`nstcalcenergy`.
2549 .. mdp:: dhdl-derivatives
2551    (yes)
2553    If yes (the default), the derivatives of the Hamiltonian with
2554    respect to lambda at each :mdp:`nstdhdl` step are written
2555    out. These values are needed for interpolation of linear energy
2556    differences with :ref:`gmx bar` (although the same can also be
2557    achieved with the right foreign lambda setting, that may not be as
2558    flexible), or with thermodynamic integration
2560 .. mdp:: dhdl-print-energy
2562    (no)
2564    Include either the total or the potential energy in the dhdl
2565    file. Options are 'no', 'potential', or 'total'. This information
2566    is needed for later free energy analysis if the states of interest
2567    are at different temperatures. If all states are at the same
2568    temperature, this information is not needed. 'potential' is useful
2569    in case one is using ``mdrun -rerun`` to generate the ``dhdl.xvg``
2570    file. When rerunning from an existing trajectory, the kinetic
2571    energy will often not be correct, and thus one must compute the
2572    residual free energy from the potential alone, with the kinetic
2573    energy component computed analytically.
2575 .. mdp:: separate-dhdl-file
2577    .. mdp-value:: yes
2579       The free energy values that are calculated (as specified with
2580       the foreign lambda and :mdp:`dhdl-derivatives` settings) are
2581       written out to a separate file, with the default name
2582       ``dhdl.xvg``. This file can be used directly with :ref:`gmx
2583       bar`.
2585    .. mdp-value:: no
2587       The free energy values are written out to the energy output file
2588       (``ener.edr``, in accumulated blocks at every :mdp:`nstenergy`
2589       steps), where they can be extracted with :ref:`gmx energy` or
2590       used directly with :ref:`gmx bar`.
2592 .. mdp:: dh-hist-size
2594    (0)
2595    If nonzero, specifies the size of the histogram into which the
2596    Delta H values (specified with foreign lambda) and the derivative
2597    dH/dl values are binned, and written to ener.edr. This can be used
2598    to save disk space while calculating free energy differences. One
2599    histogram gets written for each foreign lambda and two for the
2600    dH/dl, at every :mdp:`nstenergy` step. Be aware that incorrect
2601    histogram settings (too small size or too wide bins) can introduce
2602    errors. Do not use histograms unless you're certain you need it.
2604 .. mdp:: dh-hist-spacing
2606    (0.1)
2607    Specifies the bin width of the histograms, in energy units. Used in
2608    conjunction with :mdp:`dh-hist-size`. This size limits the
2609    accuracy with which free energies can be calculated. Do not use
2610    histograms unless you're certain you need it.
2613 Expanded Ensemble calculations
2614 ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
2616 .. mdp:: nstexpanded
2618    The number of integration steps beween attempted moves changing the
2619    system Hamiltonian in expanded ensemble simulations. Must be a
2620    multiple of :mdp:`nstcalcenergy`, but can be greater or less than
2621    :mdp:`nstdhdl`.
2623 .. mdp:: lmc-stats
2625    .. mdp-value:: no
2627       No Monte Carlo in state space is performed.
2629    .. mdp-value:: metropolis-transition
2631       Uses the Metropolis weights to update the expanded ensemble
2632       weight of each state. Min{1,exp(-(beta_new u_new - beta_old
2633       u_old)}
2635    .. mdp-value:: barker-transition
2637       Uses the Barker transition critera to update the expanded
2638       ensemble weight of each state i, defined by exp(-beta_new
2639       u_new)/(exp(-beta_new u_new)+exp(-beta_old u_old))
2641    .. mdp-value:: wang-landau
2643       Uses the Wang-Landau algorithm (in state space, not energy
2644       space) to update the expanded ensemble weights.
2646    .. mdp-value:: min-variance
2648       Uses the minimum variance updating method of Escobedo et al. to
2649       update the expanded ensemble weights. Weights will not be the
2650       free energies, but will rather emphasize states that need more
2651       sampling to give even uncertainty.
2653 .. mdp:: lmc-mc-move
2655    .. mdp-value:: no
2657       No Monte Carlo in state space is performed.
2659    .. mdp-value:: metropolis-transition
2661       Randomly chooses a new state up or down, then uses the
2662       Metropolis critera to decide whether to accept or reject:
2663       Min{1,exp(-(beta_new u_new - beta_old u_old)}
2665    .. mdp-value:: barker-transition
2667       Randomly chooses a new state up or down, then uses the Barker
2668       transition critera to decide whether to accept or reject:
2669       exp(-beta_new u_new)/(exp(-beta_new u_new)+exp(-beta_old u_old))
2671    .. mdp-value:: gibbs
2673        Uses the conditional weights of the state given the coordinate
2674        (exp(-beta_i u_i) / sum_k exp(beta_i u_i) to decide which state
2675        to move to.
2677    .. mdp-value:: metropolized-gibbs
2679        Uses the conditional weights of the state given the coordinate
2680        (exp(-beta_i u_i) / sum_k exp(beta_i u_i) to decide which state
2681        to move to, EXCLUDING the current state, then uses a rejection
2682        step to ensure detailed balance. Always more efficient that
2683        Gibbs, though only marginally so in many situations, such as
2684        when only the nearest neighbors have decent phase space
2685        overlap.
2687 .. mdp:: lmc-seed
2689    (-1)
2690    random seed to use for Monte Carlo moves in state space. When
2691    :mdp:`lmc-seed` is set to -1, a pseudo random seed is us
2693 .. mdp:: mc-temperature
2695    Temperature used for acceptance/rejection for Monte Carlo moves. If
2696    not specified, the temperature of the simulation specified in the
2697    first group of :mdp:`ref-t` is used.
2699 .. mdp:: wl-ratio
2701    (0.8)
2702    The cutoff for the histogram of state occupancies to be reset, and
2703    the free energy incrementor to be changed from delta to delta *
2704    :mdp:`wl-scale`. If we define the Nratio = (number of samples at
2705    each histogram) / (average number of samples at each
2706    histogram). :mdp:`wl-ratio` of 0.8 means that means that the
2707    histogram is only considered flat if all Nratio > 0.8 AND
2708    simultaneously all 1/Nratio > 0.8.
2710 .. mdp:: wl-scale
2712    (0.8)
2713    Each time the histogram is considered flat, then the current value
2714    of the Wang-Landau incrementor for the free energies is multiplied
2715    by :mdp:`wl-scale`. Value must be between 0 and 1.
2717 .. mdp:: init-wl-delta
2719    (1.0)
2720    The initial value of the Wang-Landau incrementor in kT. Some value
2721    near 1 kT is usually most efficient, though sometimes a value of
2722    2-3 in units of kT works better if the free energy differences are
2723    large.
2725 .. mdp:: wl-oneovert
2727    (no)
2728    Set Wang-Landau incrementor to scale with 1/(simulation time) in
2729    the large sample limit. There is significant evidence that the
2730    standard Wang-Landau algorithms in state space presented here
2731    result in free energies getting 'burned in' to incorrect values
2732    that depend on the initial state. when :mdp:`wl-oneovert` is true,
2733    then when the incrementor becomes less than 1/N, where N is the
2734    mumber of samples collected (and thus proportional to the data
2735    collection time, hence '1 over t'), then the Wang-Lambda
2736    incrementor is set to 1/N, decreasing every step. Once this occurs,
2737    :mdp:`wl-ratio` is ignored, but the weights will still stop
2738    updating when the equilibration criteria set in
2739    :mdp:`lmc-weights-equil` is achieved.
2741 .. mdp:: lmc-repeats
2743    (1)
2744    Controls the number of times that each Monte Carlo swap type is
2745    performed each iteration. In the limit of large numbers of Monte
2746    Carlo repeats, then all methods converge to Gibbs sampling. The
2747    value will generally not need to be different from 1.
2749 .. mdp:: lmc-gibbsdelta
2751    (-1)
2752    Limit Gibbs sampling to selected numbers of neighboring states. For
2753    Gibbs sampling, it is sometimes inefficient to perform Gibbs
2754    sampling over all of the states that are defined. A positive value
2755    of :mdp:`lmc-gibbsdelta` means that only states plus or minus
2756    :mdp:`lmc-gibbsdelta` are considered in exchanges up and down. A
2757    value of -1 means that all states are considered. For less than 100
2758    states, it is probably not that expensive to include all states.
2760 .. mdp:: lmc-forced-nstart
2762    (0)
2763    Force initial state space sampling to generate weights. In order to
2764    come up with reasonable initial weights, this setting allows the
2765    simulation to drive from the initial to the final lambda state,
2766    with :mdp:`lmc-forced-nstart` steps at each state before moving on
2767    to the next lambda state. If :mdp:`lmc-forced-nstart` is
2768    sufficiently long (thousands of steps, perhaps), then the weights
2769    will be close to correct. However, in most cases, it is probably
2770    better to simply run the standard weight equilibration algorithms.
2772 .. mdp:: nst-transition-matrix
2774    (-1)
2775    Frequency of outputting the expanded ensemble transition matrix. A
2776    negative number means it will only be printed at the end of the
2777    simulation.
2779 .. mdp:: symmetrized-transition-matrix
2781    (no)
2782    Whether to symmetrize the empirical transition matrix. In the
2783    infinite limit the matrix will be symmetric, but will diverge with
2784    statistical noise for short timescales. Forced symmetrization, by
2785    using the matrix T_sym = 1/2 (T + transpose(T)), removes problems
2786    like the existence of (small magnitude) negative eigenvalues.
2788 .. mdp:: mininum-var-min
2790    (100)
2791    The min-variance strategy (option of :mdp:`lmc-stats` is only
2792    valid for larger number of samples, and can get stuck if too few
2793    samples are used at each state. :mdp:`mininum-var-min` is the
2794    minimum number of samples that each state that are allowed before
2795    the min-variance strategy is activated if selected.
2797 .. mdp:: init-lambda-weights
2799    The initial weights (free energies) used for the expanded ensemble
2800    states. Default is a vector of zero weights. format is similar to
2801    the lambda vector settings in :mdp:`fep-lambdas`, except the
2802    weights can be any floating point number. Units are kT. Its length
2803    must match the lambda vector lengths.
2805 .. mdp:: lmc-weights-equil
2807    .. mdp-value:: no
2809       Expanded ensemble weights continue to be updated throughout the
2810       simulation.
2812    .. mdp-value:: yes
2814       The input expanded ensemble weights are treated as equilibrated,
2815       and are not updated throughout the simulation.
2817    .. mdp-value:: wl-delta
2819       Expanded ensemble weight updating is stopped when the
2820       Wang-Landau incrementor falls below this value.
2822    .. mdp-value:: number-all-lambda
2824       Expanded ensemble weight updating is stopped when the number of
2825       samples at all of the lambda states is greater than this value.
2827    .. mdp-value:: number-steps
2829       Expanded ensemble weight updating is stopped when the number of
2830       steps is greater than the level specified by this value.
2832    .. mdp-value:: number-samples
2834       Expanded ensemble weight updating is stopped when the number of
2835       total samples across all lambda states is greater than the level
2836       specified by this value.
2838    .. mdp-value:: count-ratio
2840       Expanded ensemble weight updating is stopped when the ratio of
2841       samples at the least sampled lambda state and most sampled
2842       lambda state greater than this value.
2844 .. mdp:: simulated-tempering
2846    (no)
2847    Turn simulated tempering on or off. Simulated tempering is
2848    implemented as expanded ensemble sampling with different
2849    temperatures instead of different Hamiltonians.
2851 .. mdp:: sim-temp-low
2853    (300) [K]
2854    Low temperature for simulated tempering.
2856 .. mdp:: sim-temp-high
2858    (300) [K]
2859    High temperature for simulated tempering.
2861 .. mdp:: simulated-tempering-scaling
2863    Controls the way that the temperatures at intermediate lambdas are
2864    calculated from the :mdp:`temperature-lambdas` part of the lambda
2865    vector.
2867    .. mdp-value:: linear
2869       Linearly interpolates the temperatures using the values of
2870       :mdp:`temperature-lambdas`, *i.e.* if :mdp:`sim-temp-low`
2871       =300, :mdp:`sim-temp-high` =400, then lambda=0.5 correspond to
2872       a temperature of 350. A nonlinear set of temperatures can always
2873       be implemented with uneven spacing in lambda.
2875    .. mdp-value:: geometric
2877       Interpolates temperatures geometrically between
2878       :mdp:`sim-temp-low` and :mdp:`sim-temp-high`. The i:th state
2879       has temperature :mdp:`sim-temp-low` * (:mdp:`sim-temp-high` /
2880       :mdp:`sim-temp-low`) raised to the power of
2881       (i/(ntemps-1)). This should give roughly equal exchange for
2882       constant heat capacity, though of course things simulations that
2883       involve protein folding have very high heat capacity peaks.
2885    .. mdp-value:: exponential
2887       Interpolates temperatures exponentially between
2888       :mdp:`sim-temp-low` and :mdp:`sim-temp-high`. The i:th state
2889       has temperature :mdp:`sim-temp-low` + (:mdp:`sim-temp-high` -
2890       :mdp:`sim-temp-low`)*((exp(:mdp:`temperature-lambdas`
2891       (i))-1)/(exp(1.0)-i)).
2894 Non-equilibrium MD
2895 ^^^^^^^^^^^^^^^^^^
2897 .. mdp:: acc-grps
2899    groups for constant acceleration (*e.g.* ``Protein Sol``) all atoms
2900    in groups Protein and Sol will experience constant acceleration as
2901    specified in the :mdp:`accelerate` line
2903 .. mdp:: accelerate
2905    (0) [nm ps\ :sup:`-2`]
2906    acceleration for :mdp:`acc-grps`; x, y and z for each group
2907    (*e.g.* ``0.1 0.0 0.0 -0.1 0.0 0.0`` means that first group has
2908    constant acceleration of 0.1 nm ps\ :sup:`-2` in X direction, second group
2909    the opposite).
2911 .. mdp:: freezegrps
2913    Groups that are to be frozen (*i.e.* their X, Y, and/or Z position
2914    will not be updated; *e.g.* ``Lipid SOL``). :mdp:`freezedim`
2915    specifies for which dimension(s) the freezing applies. To avoid
2916    spurious contributions to the virial and pressure due to large
2917    forces between completely frozen atoms you need to use energy group
2918    exclusions, this also saves computing time. Note that coordinates
2919    of frozen atoms are not scaled by pressure-coupling algorithms.
2921 .. mdp:: freezedim
2923    dimensions for which groups in :mdp:`freezegrps` should be frozen,
2924    specify ``Y`` or ``N`` for X, Y and Z and for each group (*e.g.*
2925    ``Y Y N N N N`` means that particles in the first group can move only in
2926    Z direction. The particles in the second group can move in any
2927    direction).
2929 .. mdp:: cos-acceleration
2931    (0) [nm ps\ :sup:`-2`]
2932    the amplitude of the acceleration profile for calculating the
2933    viscosity. The acceleration is in the X-direction and the magnitude
2934    is :mdp:`cos-acceleration` cos(2 pi z/boxheight). Two terms are
2935    added to the energy file: the amplitude of the velocity profile and
2936    1/viscosity.
2938 .. mdp:: deform
2940    (0 0 0 0 0 0) [nm ps\ :sup:`-1`]
2941    The velocities of deformation for the box elements: a(x) b(y) c(z)
2942    b(x) c(x) c(y). Each step the box elements for which :mdp:`deform`
2943    is non-zero are calculated as: box(ts)+(t-ts)*deform, off-diagonal
2944    elements are corrected for periodicity. The coordinates are
2945    transformed accordingly. Frozen degrees of freedom are (purposely)
2946    also transformed. The time ts is set to t at the first step and at
2947    steps at which x and v are written to trajectory to ensure exact
2948    restarts. Deformation can be used together with semiisotropic or
2949    anisotropic pressure coupling when the appropriate
2950    compressibilities are set to zero. The diagonal elements can be
2951    used to strain a solid. The off-diagonal elements can be used to
2952    shear a solid or a liquid.
2955 Electric fields
2956 ^^^^^^^^^^^^^^^
2958 .. mdp:: electric-field-x
2959 .. mdp:: electric-field-y
2960 .. mdp:: electric-field-z
2962    Here you can specify an electric field that optionally can be
2963    alternating and pulsed. The general expression for the field
2964    has the form of a gaussian laser pulse:
2966    .. math:: E(t) = E_0 \exp\left[-\frac{(t-t_0)^2}{2\sigma^2}\right]\cos\left[\omega (t-t_0)\right]
2968    For example, the four parameters for direction x are set in the
2969    fields of :mdp:`electric-field-x` (and similar for ``electric-field-y``
2970    and ``electric-field-z``) like
2972    ``electric-field-x  = E0 omega t0 sigma``
2974    with units (respectively) V nm\ :sup:`-1`, ps\ :sup:`-1`, ps, ps.
2976    In the special case that ``sigma = 0``, the exponential term is omitted
2977    and only the cosine term is used. If also ``omega = 0`` a static
2978    electric field is applied.
2980    Read more at :ref:`electric fields` and in ref. \ :ref:`146 <refCaleman2008a>`.
2983 Mixed quantum/classical molecular dynamics
2984 ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
2986 .. mdp:: QMMM-grps
2988    groups to be descibed at the QM level for MiMiC QM/MM
2990 .. MDP:: QMMM
2992    .. mdp-value:: no
2994       QM/MM is no longer supported via these .mdp options. For MiMic, use no here.
2996 Computational Electrophysiology
2997 ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
2998 Use these options to switch on and control ion/water position exchanges in "Computational
2999 Electrophysiology" simulation setups. (See the `reference manual`_ for details).
3001 .. mdp:: swapcoords
3003    .. mdp-value:: no
3005       Do not enable ion/water position exchanges.
3007    .. mdp-value:: X ; Y ; Z
3009       Allow for ion/water position exchanges along the chosen direction.
3010       In a typical setup with the membranes parallel to the x-y plane,
3011       ion/water pairs need to be exchanged in Z direction to sustain the
3012       requested ion concentrations in the compartments.
3014 .. mdp:: swap-frequency
3016    (1) The swap attempt frequency, i.e. every how many time steps the ion counts
3017    per compartment are determined and exchanges made if necessary.
3018    Normally it is not necessary to check at every time step.
3019    For typical Computational Electrophysiology setups, a value of about 100 is
3020    sufficient and yields a negligible performance impact.
3022 .. mdp:: split-group0
3024    Name of the index group of the membrane-embedded part of channel #0.
3025    The center of mass of these atoms defines one of the compartment boundaries
3026    and should be chosen such that it is near the center of the membrane.
3028 .. mdp:: split-group1
3030    Channel #1 defines the position of the other compartment boundary.
3032 .. mdp:: massw-split0
3034    (no) Defines whether or not mass-weighting is used to calculate the split group center.
3036    .. mdp-value:: no
3038       Use the geometrical center.
3040    .. mdp-value:: yes
3042       Use the center of mass.
3044 .. mdp:: massw-split1
3046    (no) As above, but for split-group #1.
3048 .. mdp:: solvent-group
3050    Name of the index group of solvent molecules.
3052 .. mdp:: coupl-steps
3054    (10) Average the number of ions per compartment over these many swap attempt steps.
3055    This can be used to prevent that ions near a compartment boundary
3056    (diffusing through a channel, e.g.) lead to unwanted back and forth swaps.
3058 .. mdp:: iontypes
3060    (1) The number of different ion types to be controlled. These are during the
3061    simulation exchanged with solvent molecules to reach the desired reference numbers.
3063 .. mdp:: iontype0-name
3065    Name of the first ion type.
3067 .. mdp:: iontype0-in-A
3069    (-1) Requested (=reference) number of ions of type 0 in compartment A.
3070    The default value of -1 means: use the number of ions as found in time step 0
3071    as reference value.
3073 .. mdp:: iontype0-in-B
3075    (-1) Reference number of ions of type 0 for compartment B.
3077 .. mdp:: bulk-offsetA
3079    (0.0) Offset of the first swap layer from the compartment A midplane.
3080    By default (i.e. bulk offset = 0.0), ion/water exchanges happen between layers
3081    at maximum distance (= bulk concentration) to the split group layers. However,
3082    an offset b (-1.0 < b < +1.0) can be specified to offset the bulk layer from the middle at 0.0
3083    towards one of the compartment-partitioning layers (at +/- 1.0).
3085 .. mdp:: bulk-offsetB
3087    (0.0) Offset of the other swap layer from the compartment B midplane.
3090 .. mdp:: threshold
3092    (\1) Only swap ions if threshold difference to requested count is reached.
3094 .. mdp:: cyl0-r
3096    (2.0) [nm] Radius of the split cylinder #0.
3097    Two split cylinders (mimicking the channel pores) can optionally be defined
3098    relative to the center of the split group. With the help of these cylinders
3099    it can be counted which ions have passed which channel. The split cylinder
3100    definition has no impact on whether or not ion/water swaps are done.
3102 .. mdp:: cyl0-up
3104    (1.0) [nm] Upper extension of the split cylinder #0.
3106 .. mdp:: cyl0-down
3108    (1.0) [nm] Lower extension of the split cylinder #0.
3110 .. mdp:: cyl1-r
3112    (2.0) [nm] Radius of the split cylinder #1.
3114 .. mdp:: cyl1-up
3116    (1.0) [nm] Upper extension of the split cylinder #1.
3118 .. mdp:: cyl1-down
3120    (1.0) [nm] Lower extension of the split cylinder #1.
3122 Density-guided simulations
3123 ^^^^^^^^^^^^^^^^^^^^^^^^^^
3125 These options enable and control the calculation and application of additional
3126 forces that are derived from three-dimensional densities, e.g., from cryo
3127 electron-microscopy experiments. (See the `reference manual`_ for details)
3129 .. mdp:: density-guided-simulation-active
3131    (no) Activate density-guided simulations.
3133 .. mdp:: density-guided-simulation-group
3135    (protein) The atoms that are subject to the forces from the density-guided
3136    simulation and contribute to the simulated density.
3138 .. mdp:: density-guided-simulation-similarity-measure
3140    (inner-product) Similarity measure between the density that is calculated
3141    from the atom positions and the reference density.
3143    .. mdp-value:: inner-product
3145       Takes the sum of the product of reference density and simulated density
3146       voxel values.
3148    .. mdp-value:: relative-entropy
3150       Uses the negative relative entropy (or Kullback-Leibler divergence)
3151       between reference density and simulated density as similarity measure.
3152       Negative density values are ignored.
3154    .. mdp-value:: cross-correlation
3156       Uses the Pearson correlation coefficient between reference density and
3157       simulated density as similarity measure.
3159 .. mdp:: density-guided-simulation-atom-spreading-weight
3161    (unity) Determines the multiplication factor for the Gaussian kernel when
3162    spreading atoms on the grid.
3164    .. mdp-value:: unity
3166       Every atom in the density fitting group is assigned the same unit factor.
3168    .. mdp-value:: mass
3170       Atoms contribute to the simulated density proportional to their mass.
3172    .. mdp-value:: charge
3174       Atoms contribute to the simulated density proportional to their charge.
3176 .. mdp:: density-guided-simulation-force-constant
3178    (1e+09) [kJ mol\ :sup:`-1`] The scaling factor for density-guided simulation
3179    forces. May also be negative.
3181 .. mdp:: density-guided-simulation-gaussian-transform-spreading-width
3183    (0.2) [nm] The Gaussian RMS width for the spread kernel for the simulated
3184    density.
3186 .. mdp:: density-guided-simulation-gaussian-transform-spreading-range-in-multiples-of-width
3188    (4) The range after which the gaussian is cut off in multiples of the Gaussian
3189    RMS width described above.
3191 .. mdp:: density-guided-simulation-reference-density-filename
3193    (reference.mrc) Reference density file name using an absolute path or a path
3194    relative to the to the folder from which :ref:`gmx mdrun` is called.
3196 .. mdp:: density-guided-simulation-nst
3198    (1) Interval in steps at which the density fitting forces are evaluated
3199    and applied. The forces are scaled by this number when applied (See the
3200    `reference manual`_ for details).
3202 .. mdp:: density-guided-simulation-normalize-densities
3204    (true) Normalize the sum of density voxel values to one for the reference
3205    density as well as the simulated density.
3207 .. mdp:: density-guided-simulation-adaptive-force-scaling
3209    (false) Adapt the force constant to ensure a steady increase in similarity
3210    between simulated and reference density.
3212    .. mdp-value: false
3214       Do not use adaptive force scaling.
3216    .. mdp-value:: true
3218       Use adaptive force scaling.
3220 .. mdp:: density-guided-simulation-adaptive-force-scaling-time-constant
3222    (4) [ps] Couple force constant to increase in similarity with reference density
3223    with this time constant. Larger times result in looser coupling.
3225 User defined thingies
3226 ^^^^^^^^^^^^^^^^^^^^^
3228 .. mdp:: user1-grps
3229 .. mdp:: user2-grps
3230 .. mdp:: userint1 (0)
3231 .. mdp:: userint2 (0)
3232 .. mdp:: userint3 (0)
3233 .. mdp:: userint4 (0)
3234 .. mdp:: userreal1 (0)
3235 .. mdp:: userreal2 (0)
3236 .. mdp:: userreal3 (0)
3237 .. mdp:: userreal4 (0)
3239    These you can use if you modify code. You can pass integers and
3240    reals and groups to your subroutine. Check the inputrec definition
3241    in ``src/gromacs/mdtypes/inputrec.h``
3243 Removed features
3244 ^^^^^^^^^^^^^^^^
3246 These features have been removed from |Gromacs|, but so that old
3247 :ref:`mdp` and :ref:`tpr` files cannot be mistakenly misused, we still
3248 parse this option. :ref:`gmx grompp` and :ref:`gmx mdrun` will issue a
3249 fatal error if this is set.
3251 .. mdp:: adress
3253    (no)
3255 .. mdp:: implicit-solvent
3257    (no)
3259 .. _reference manual: gmx-manual-parent-dir_