Add membrane embedding docs
[alexxy/gromacs.git] / docs / user-guide / mdrun-features.rst
1 Useful mdrun features
2 =======================
3 This section discusses features in :ref:`gmx mdrun` that don't fit well
4 elsewhere.
5
6 .. _single-point energy:
7
8 Re-running a simulation
9 -----------------------
10 The rerun feature allows you to take any trajectory file ``traj.trr``
11 and compute quantities based upon the coordinates in that file using
12 the model physics supplied in the ``topol.tpr`` file. It can be used
13 with command lines like ``mdrun -s topol -rerun traj.trr``. That :ref:`tpr`
14 could be different from the one that generated the trajectory. This
15 can be used to compute the energy or forces for exactly the
16 coordinates supplied as input, or to extract quantities based on
17 subsets of the molecular system (see :ref:`gmx convert-tpr` and
18 :ref:`gmx trjconv`). It is easier to do a correct "single-point" energy
19 evaluation with this feature than a 0-step simulation.
20
21 Neighbour searching is normally performed for every frame in the
22 trajectory, since :ref:`gmx mdrun` can no longer assume anything about how the
23 structures were generated. If :mdp:`nstlist` is zero, then only one
24 neighbour list will be constructed. Naturally, no update or constraint
25 algorithms are ever used.
26
27 Running a simulation in reproducible mode
28 -----------------------------------------
29 It is generally difficult to run an efficient parallel MD simulation
30 that is based primarily on floating-point arithmetic and is fully
31 reproducible. By default, :ref:`gmx mdrun` will observe how things are going
32 and vary how the simulation is conducted in order to optimize
33 throughput. However, there is a "reproducible mode" available with
34 ``mdrun -reprod`` that will systematically eliminate all sources of
35 variation within that run; repeated invocations on the same input and
36 hardware will be binary identical. However, running in this mode on
37 different hardware, or with a different compiler, etc. will not be
38 reproducible. This should normally only be used when investigating
39 possible problems.
40
41 Running multi-simulations
42 -------------------------
43 There are numerous situations where running a related set of
44 simulations within the same invocation of mdrun are necessary or
45 useful. Running a replica-exchange simulation requires it, as do
46 simulations using ensemble-based distance or orientation restraints.
47 Running a related series of lambda points for a free-energy
48 computation is also convenient to do this way.
49
50 This feature requires
51 :ref:`configuring GROMACS with an external MPI library <mpi-support>`
52 so that the set of
53 simulations can communicate. The ``n`` simulations within the set can
54 use internal MPI parallelism also, so that ``mpirun -np x mdrun_mpi``
55 for ``x`` a multiple of ``n`` will use ``x/n`` ranks per simulation.
56
57 There are two ways of organizing files when running such
58 simulations. All of the normal mechanisms work in either case,
59 including ``-deffnm``.
60
61 ``-multidir``
62    You must create a set of ``n`` directories for the ``n`` simulations,
63    place all the relevant input files in those directories (e.g. named
64    ``topol.tpr``), and run with
65    ``mpirun -np x mdrun_mpi -s topol -multidir <names-of-directories>``.
66    If the order of the simulations
67    within the multi-simulation is significant, then you are responsible
68    for ordering their names when you provide them to ``-multidir``. Be
69    careful with shells that do filename globbing dictionary-style, e.g.
70    ``dir1 dir10 dir11 ... dir2 ...``. This option is generally the
71    most convenient to use. ``mdrun -table`` for the group cutoff-scheme
72    works only in this mode.
73
74 ``-multi``
75    You must organize that the filenames for each simulation in a set of
76    ``n`` simulations have an integer ``0`` through ``n-1`` appended to
77    the filename (e.g. ``topol2.tpr``), and run with
78    ``mpirun -np x mdrun -multi n -s input``. The order of simulations
79    within the set is determined by the integer appended.
80
81 Examples running multi-simulations
82 ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
83
84 ::
85
86     mpirun -np 32 mdrun_mpi -multi
87
88 Starts a multi-simulation on 32 ranks with as many simulations ``n`` as
89 there are files named ``topol*.tpr`` for integers ``0`` to ``n-1``. Other
90 input and output files are suffixed similarly.
91
92 ::
93
94     mpirun -np 32 mdrun_mpi -multidir a b c d
95
96 Starts a multi-simulation on 32 ranks with 4 simulations. The input
97 and output files are found in directories ``a``, ``b``, ``c``, and ``d``.
98
99 ::
100
101     mpirun -np 32 mdrun_mpi -multidir a b c d -gpu_id 0000000011111111
102
103 Starts the same multi-simulation as before. On a machine with two
104 physical nodes and two GPUs per node, there will be 16 MPI ranks per
105 node, and 8 MPI ranks per simulation. The 16 MPI ranks doing PP work
106 on a node are mapped to the GPUs with IDs 0 and 1, even though they
107 come from more than one simulation. They are mapped in the order
108 indicated, so that the PP ranks from each simulation use a single
109 GPU. However, the order ``0101010101010101`` could run faster.
110
111 Running replica-exchange simulations
112 ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
113
114 When running a multi-simulation, using ``mdrun -replex n`` means that a
115 replica exchange is attempted every given number of steps. The number
116 of replicas is set with the ``-multi`` or ``-multidir`` option, described
117 above.  All run input files should use a different value for the
118 coupling parameter (e.g. temperature), which ascends over the set of
119 input files. The random seed for replica exchange is set with
120 ``-reseed``. After every exchange, the velocities are scaled and
121 neighbor searching is performed. See the Reference Manual for more
122 details on how replica exchange functions in GROMACS.
123
124 Controlling the length of the simulation
125 ----------------------------------------
126
127 Normally, the length of an MD simulation is best managed through the
128 [.mdp] option [nsteps](#nsteps), however there are situations where
129 more control is useful. `mdrun -nsteps 100` overrides the [.mdp] file
130 and executes 100 steps. `mdrun -maxh 2.5` will terminate the
131 simulation shortly before 2.5 hours elapse, which can be useful when
132 running under cluster queues (as long as the queueing system does not
133 ever suspend the simulation).
134
135 Running a membrane protein embedding simulation
136 -----------------------------------------------
137
138 This is a module to help embed a membrane protein into an equilibrated
139 lipid bilayer at a position and orientation specified by the user. The
140 main advantage is that it is possible to use very complex lipid bilayers
141 with a number of different components that have been relaxed for a
142 long time in a previous simulation. In theory that could be accomplished
143 with a procedure similar to genbox, but since lipids are much larger
144 than water molecules it will lead to a large vacuum layer between the
145 protein and membrane if we remove all molecules where any atom is
146 overlapping. Instead, this module works by first artifically shrinking
147 the protein in the xy-plane, then it removes lipids that overlap with
148 a much smaller core, after which we gradually push the protein atoms
149 back to their initial positions, while using normal dynamics for the
150 rest of the system so lipids adapt to the protein.
151
152 To use membrane embedding, start by building a lipid bilayer that is
153 just-so-slightly larger in the xy-plane than what you expect to need
154 in the end, and make sure you have enough water outside the membrane
155 to accomodate globular domains. Place the protein in the same coordinate
156 file (and topology) as your lipid bilayer, and make sure it is in the
157 orientation and position you want right in the middle of the bilayer.
158
159 The first settings have to be entered in the mdp file that controls
160 your simulation. You need an energy group corresponding to your
161 protein, this group should be frozen (all dimensions), and we should
162 exclude all interactions inside the protein to avoid problems when it
163 is distorted. For instance:
164
165 ::
166
167     integrator     = md
168     energygrps     = Protein
169     freezegrps     = Protein
170     freezedim      = Y Y Y
171     energygrp_excl = Protein Protein
172
173 You will also need a number of settings for the actual membrane
174 embedding process. These are entered as similar name and value pairs,
175 but in the separate text data file ``embed.dat`` that you provide as
176 the argument to the ``-membed`` option (we refer to the below
177 when explaining the process). The embedding works in for stages:
178
179 1. The protein is resized around its center of mass by a factor
180    ``xy`` in the xy-plane (the bilayer plane), and a factor ``z``
181    along the z-axis (normal to the bilayer). If the height of the
182    protein is the same or smaller than the thickness of the
183    membrane, a z-fraction larger than 1.0 can prevent the protein
184    from being enveloped by the lipids.
185
186 2. All lipid and solvent molecules overlapping with the resized
187    protein are removed. All interactions inside the protein are
188    turned off to prevent numerical issues for small values of the
189    scaling fractions.
190
191 3. A single md step is performed, where atoms in the rest of the
192    system are moved.
193
194 4. The resize factors are adjusted by the small amounts
195    (1-xy)/nxy and (1-z)/nz, where ``nxy`` and ``nz`` are the
196    number of iterations to use.  The resize factor for the xy-plane
197    is adjusted first. The resize factor for the z-direction is not
198    changed until the xy factor is 1.0 (after ``nxy`` iterations).
199
200 5. Steps 3 and 4 are repeated until the protein has again reached
201    its original size, i.e. after nxy+nz iterations. After the
202    embedding you might still want to perform a short relaxation.
203
204 Parameters that can be specified in ``embed.dat``, with default
205 values that will be used if the setting is omitted:
206
207 - ``xyinit`` (0.5) Resize factor for the protein in the xy
208   dimension before starting embedding.
209
210 - ``xyend`` (1.0) Final resize factor in the xy dimension.
211
212 - ``zinit`` (1.0) Resize factor for the protein in the z
213    dimension before starting embedding.
214
215 - ``zend`` (1.0) Final resize faction in the z dimension.
216
217 - ``nxy`` (1000) Number of iteration for the xy dimension.
218
219 - ``nz`` (0) Number of iterations for the z dimension.
220
221 - ``rad`` (0.22) Probe radius to check for overlap between
222   the group to embed and the membrane.
223
224 - ``pieces`` (1) Perform piecewise resize. Select parts of the group
225   to insert and resize these with respect to their own geometrical center.
226
227 - ``asymmetry`` (no) Allow asymmetric insertion, i.e. the number of
228   lipids removed from the upper and lower leaflet will not be checked.
229
230 - ``ndiff`` (0) Number of lipids that will additionally be removed
231   from the lower (negative number) or upper (positive number)
232   membrane leaflet.
233
234 - ``maxwarn`` (0) Largest number of membed warnings allowed.