Clarified docs for use of mdrun -nstlist
[alexxy/gromacs.git] / src / programs / mdrun / mdrun.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2011,2012,2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "mdrun_main.h"
38
39 #ifdef HAVE_CONFIG_H
40 #include "config.h"
41 #endif
42
43 #include <stdio.h>
44
45 #include "gromacs/legacyheaders/checkpoint.h"
46 #include "gromacs/legacyheaders/copyrite.h"
47 #include "gromacs/legacyheaders/gmx_fatal.h"
48 #include "gromacs/legacyheaders/macros.h"
49 #include "gromacs/legacyheaders/main.h"
50 #include "gromacs/legacyheaders/mdrun.h"
51 #include "gromacs/legacyheaders/network.h"
52 #include "gromacs/legacyheaders/readinp.h"
53 #include "gromacs/legacyheaders/typedefs.h"
54 #include "gromacs/legacyheaders/types/commrec.h"
55
56 #include "gromacs/commandline/pargs.h"
57 #include "gromacs/fileio/filenm.h"
58
59 int gmx_mdrun(int argc, char *argv[])
60 {
61     const char   *desc[] = {
62         "[THISMODULE] is the main computational chemistry engine",
63         "within GROMACS. Obviously, it performs Molecular Dynamics simulations,",
64         "but it can also perform Stochastic Dynamics, Energy Minimization,",
65         "test particle insertion or (re)calculation of energies.",
66         "Normal mode analysis is another option. In this case [TT]mdrun[tt]",
67         "builds a Hessian matrix from single conformation.",
68         "For usual Normal Modes-like calculations, make sure that",
69         "the structure provided is properly energy-minimized.",
70         "The generated matrix can be diagonalized by [gmx-nmeig].[PAR]",
71         "The [TT]mdrun[tt] program reads the run input file ([TT]-s[tt])",
72         "and distributes the topology over nodes if needed.",
73         "[TT]mdrun[tt] produces at least four output files.",
74         "A single log file ([TT]-g[tt]) is written, unless the option",
75         "[TT]-seppot[tt] is used, in which case each node writes a log file.",
76         "The trajectory file ([TT]-o[tt]), contains coordinates, velocities and",
77         "optionally forces.",
78         "The structure file ([TT]-c[tt]) contains the coordinates and",
79         "velocities of the last step.",
80         "The energy file ([TT]-e[tt]) contains energies, the temperature,",
81         "pressure, etc, a lot of these things are also printed in the log file.",
82         "Optionally coordinates can be written to a compressed trajectory file",
83         "([TT]-x[tt]).[PAR]",
84         "The option [TT]-dhdl[tt] is only used when free energy calculation is",
85         "turned on.[PAR]",
86         "A simulation can be run in parallel using two different parallelization",
87         "schemes: MPI parallelization and/or OpenMP thread parallelization.",
88         "The MPI parallelization uses multiple processes when [TT]mdrun[tt] is",
89         "compiled with a normal MPI library or threads when [TT]mdrun[tt] is",
90         "compiled with the GROMACS built-in thread-MPI library. OpenMP threads",
91         "are supported when [TT]mdrun[tt] is compiled with OpenMP. Full OpenMP support",
92         "is only available with the Verlet cut-off scheme, with the (older)",
93         "group scheme only PME-only processes can use OpenMP parallelization.",
94         "In all cases [TT]mdrun[tt] will by default try to use all the available",
95         "hardware resources. With a normal MPI library only the options",
96         "[TT]-ntomp[tt] (with the Verlet cut-off scheme) and [TT]-ntomp_pme[tt],",
97         "for PME-only processes, can be used to control the number of threads.",
98         "With thread-MPI there are additional options [TT]-nt[tt], which sets",
99         "the total number of threads, and [TT]-ntmpi[tt], which sets the number",
100         "of thread-MPI threads.",
101         "The number of OpenMP threads used by [TT]mdrun[tt] can also be set with",
102         "the standard environment variable, [TT]OMP_NUM_THREADS[tt].",
103         "The [TT]GMX_PME_NUM_THREADS[tt] environment variable can be used to specify",
104         "the number of threads used by the PME-only processes.[PAR]",
105         "Note that combined MPI+OpenMP parallelization is in many cases",
106         "slower than either on its own. However, at high parallelization, using the",
107         "combination is often beneficial as it reduces the number of domains and/or",
108         "the number of MPI ranks. (Less and larger domains can improve scaling,",
109         "with separate PME processes fewer MPI ranks reduces communication cost.)",
110         "OpenMP-only parallelization is typically faster than MPI-only parallelization",
111         "on a single CPU(-die). Since we currently don't have proper hardware",
112         "topology detection, [TT]mdrun[tt] compiled with thread-MPI will only",
113         "automatically use OpenMP-only parallelization when you use up to 4",
114         "threads, up to 12 threads with Intel Nehalem/Westmere, or up to 16",
115         "threads with Intel Sandy Bridge or newer CPUs. Otherwise MPI-only",
116         "parallelization is used (except with GPUs, see below).",
117         "[PAR]",
118         "To quickly test the performance of the new Verlet cut-off scheme",
119         "with old [TT].tpr[tt] files, either on CPUs or CPUs+GPUs, you can use",
120         "the [TT]-testverlet[tt] option. This should not be used for production,",
121         "since it can slightly modify potentials and it will remove charge groups",
122         "making analysis difficult, as the [TT].tpr[tt] file will still contain",
123         "charge groups. For production simulations it is highly recommended",
124         "to specify [TT]cutoff-scheme = Verlet[tt] in the [TT].mdp[tt] file.",
125         "[PAR]",
126         "With GPUs (only supported with the Verlet cut-off scheme), the number",
127         "of GPUs should match the number of MPI processes or MPI threads,",
128         "excluding PME-only processes/threads. With thread-MPI, unless set on the command line, the number",
129         "of MPI threads will automatically be set to the number of GPUs detected.",
130         "To use a subset of the available GPUs, or to manually provide a mapping of",
131         "GPUs to PP ranks, you can use the [TT]-gpu_id[tt] option. The argument of [TT]-gpu_id[tt] is",
132         "a string of digits (without delimiter) representing device id-s of the GPUs to be used.",
133         "For example, \"[TT]02[tt]\" specifies using GPUs 0 and 2 in the first and second PP ranks per compute node",
134         "respectively. To select different sets of GPU-s",
135         "on different nodes of a compute cluster, use the [TT]GMX_GPU_ID[tt] environment",
136         "variable instead. The format for [TT]GMX_GPU_ID[tt] is identical to ",
137         "[TT]-gpu_id[tt], with the difference that an environment variable can have",
138         "different values on different compute nodes. Multiple MPI ranks on each node",
139         "can share GPUs. This is accomplished by specifying the id(s) of the GPU(s)",
140         "multiple times, e.g. \"[TT]0011[tt]\" for four ranks sharing two GPUs in this node.",
141         "This works within a single simulation, or a multi-simulation, with any form of MPI.",
142         "[PAR]",
143         "With the Verlet cut-off scheme and verlet-buffer-tolerance set,",
144         "the pair-list update interval nstlist can be chosen freely with",
145         "the option [TT]-nstlist[tt]. [TT]mdrun[tt] will then adjust",
146         "the pair-list cut-off to maintain accuracy, and not adjust nstlist.",
147         "Otherwise, by default, [TT]mdrun[tt] will try to increase the",
148         "value of nstlist set in the [TT].mdp[tt] file to improve the",
149         "performance. For CPU-only runs, nstlist might increase to 20, for",
150         "GPU runs up to 40. For medium to high parallelization or with",
151         "fast GPUs, a (user-supplied) larger nstlist value can give much",
152         "better performance.",
153         "[PAR]",
154         "When using PME with separate PME nodes or with a GPU, the two major",
155         "compute tasks, the non-bonded force calculation and the PME calculation",
156         "run on different compute resources. If this load is not balanced,",
157         "some of the resources will be idle part of time. With the Verlet",
158         "cut-off scheme this load is automatically balanced when the PME load",
159         "is too high (but not when it is too low). This is done by scaling",
160         "the Coulomb cut-off and PME grid spacing by the same amount. In the first",
161         "few hundred steps different settings are tried and the fastest is chosen",
162         "for the rest of the simulation. This does not affect the accuracy of",
163         "the results, but it does affect the decomposition of the Coulomb energy",
164         "into particle and mesh contributions. The auto-tuning can be turned off",
165         "with the option [TT]-notunepme[tt].",
166         "[PAR]",
167         "[TT]mdrun[tt] pins (sets affinity of) threads to specific cores,",
168         "when all (logical) cores on a compute node are used by [TT]mdrun[tt],",
169         "even when no multi-threading is used,",
170         "as this usually results in significantly better performance.",
171         "If the queuing systems or the OpenMP library pinned threads, we honor",
172         "this and don't pin again, even though the layout may be sub-optimal.",
173         "If you want to have [TT]mdrun[tt] override an already set thread affinity",
174         "or pin threads when using less cores, use [TT]-pin on[tt].",
175         "With SMT (simultaneous multithreading), e.g. Intel Hyper-Threading,",
176         "there are multiple logical cores per physical core.",
177         "The option [TT]-pinstride[tt] sets the stride in logical cores for",
178         "pinning consecutive threads. Without SMT, 1 is usually the best choice.",
179         "With Intel Hyper-Threading 2 is best when using half or less of the",
180         "logical cores, 1 otherwise. The default value of 0 do exactly that:",
181         "it minimizes the threads per logical core, to optimize performance.",
182         "If you want to run multiple [TT]mdrun[tt] jobs on the same physical node,"
183         "you should set [TT]-pinstride[tt] to 1 when using all logical cores.",
184         "When running multiple [TT]mdrun[tt] (or other) simulations on the same physical",
185         "node, some simulations need to start pinning from a non-zero core",
186         "to avoid overloading cores; with [TT]-pinoffset[tt] you can specify",
187         "the offset in logical cores for pinning.",
188         "[PAR]",
189         "When [TT]mdrun[tt] is started using MPI with more than 1 process",
190         "or with thread-MPI with more than 1 thread, MPI parallelization is used.",
191         "Domain decomposition is always used with MPI parallelism.",
192         "[PAR]",
193         "With domain decomposition, the spatial decomposition can be set",
194         "with option [TT]-dd[tt]. By default [TT]mdrun[tt] selects a good decomposition.",
195         "The user only needs to change this when the system is very inhomogeneous.",
196         "Dynamic load balancing is set with the option [TT]-dlb[tt],",
197         "which can give a significant performance improvement,",
198         "especially for inhomogeneous systems. The only disadvantage of",
199         "dynamic load balancing is that runs are no longer binary reproducible,",
200         "but in most cases this is not important.",
201         "By default the dynamic load balancing is automatically turned on",
202         "when the measured performance loss due to load imbalance is 5% or more.",
203         "At low parallelization these are the only important options",
204         "for domain decomposition.",
205         "At high parallelization the options in the next two sections",
206         "could be important for increasing the performace.",
207         "[PAR]",
208         "When PME is used with domain decomposition, separate nodes can",
209         "be assigned to do only the PME mesh calculation;",
210         "this is computationally more efficient starting at about 12 nodes.",
211         "The number of PME nodes is set with option [TT]-npme[tt],",
212         "this can not be more than half of the nodes.",
213         "By default [TT]mdrun[tt] makes a guess for the number of PME",
214         "nodes when the number of nodes is larger than 11 or performance wise",
215         "not compatible with the PME grid x dimension.",
216         "But the user should optimize npme. Performance statistics on this issue",
217         "are written at the end of the log file.",
218         "For good load balancing at high parallelization, the PME grid x and y",
219         "dimensions should be divisible by the number of PME nodes",
220         "(the simulation will run correctly also when this is not the case).",
221         "[PAR]",
222         "This section lists all options that affect the domain decomposition.",
223         "[PAR]",
224         "Option [TT]-rdd[tt] can be used to set the required maximum distance",
225         "for inter charge-group bonded interactions.",
226         "Communication for two-body bonded interactions below the non-bonded",
227         "cut-off distance always comes for free with the non-bonded communication.",
228         "Atoms beyond the non-bonded cut-off are only communicated when they have",
229         "missing bonded interactions; this means that the extra cost is minor",
230         "and nearly indepedent of the value of [TT]-rdd[tt].",
231         "With dynamic load balancing option [TT]-rdd[tt] also sets",
232         "the lower limit for the domain decomposition cell sizes.",
233         "By default [TT]-rdd[tt] is determined by [TT]mdrun[tt] based on",
234         "the initial coordinates. The chosen value will be a balance",
235         "between interaction range and communication cost.",
236         "[PAR]",
237         "When inter charge-group bonded interactions are beyond",
238         "the bonded cut-off distance, [TT]mdrun[tt] terminates with an error message.",
239         "For pair interactions and tabulated bonds",
240         "that do not generate exclusions, this check can be turned off",
241         "with the option [TT]-noddcheck[tt].",
242         "[PAR]",
243         "When constraints are present, option [TT]-rcon[tt] influences",
244         "the cell size limit as well.",
245         "Atoms connected by NC constraints, where NC is the LINCS order plus 1,",
246         "should not be beyond the smallest cell size. A error message is",
247         "generated when this happens and the user should change the decomposition",
248         "or decrease the LINCS order and increase the number of LINCS iterations.",
249         "By default [TT]mdrun[tt] estimates the minimum cell size required for P-LINCS",
250         "in a conservative fashion. For high parallelization it can be useful",
251         "to set the distance required for P-LINCS with the option [TT]-rcon[tt].",
252         "[PAR]",
253         "The [TT]-dds[tt] option sets the minimum allowed x, y and/or z scaling",
254         "of the cells with dynamic load balancing. [TT]mdrun[tt] will ensure that",
255         "the cells can scale down by at least this factor. This option is used",
256         "for the automated spatial decomposition (when not using [TT]-dd[tt])",
257         "as well as for determining the number of grid pulses, which in turn",
258         "sets the minimum allowed cell size. Under certain circumstances",
259         "the value of [TT]-dds[tt] might need to be adjusted to account for",
260         "high or low spatial inhomogeneity of the system.",
261         "[PAR]",
262         "The option [TT]-gcom[tt] can be used to only do global communication",
263         "every n steps.",
264         "This can improve performance for highly parallel simulations",
265         "where this global communication step becomes the bottleneck.",
266         "For a global thermostat and/or barostat the temperature",
267         "and/or pressure will also only be updated every [TT]-gcom[tt] steps.",
268         "By default it is set to the minimum of nstcalcenergy and nstlist.[PAR]",
269         "With [TT]-rerun[tt] an input trajectory can be given for which ",
270         "forces and energies will be (re)calculated. Neighbor searching will be",
271         "performed for every frame, unless [TT]nstlist[tt] is zero",
272         "(see the [TT].mdp[tt] file).[PAR]",
273         "ED (essential dynamics) sampling and/or additional flooding potentials",
274         "are switched on by using the [TT]-ei[tt] flag followed by an [TT].edi[tt]",
275         "file. The [TT].edi[tt] file can be produced with the [TT]make_edi[tt] tool",
276         "or by using options in the essdyn menu of the WHAT IF program.",
277         "[TT]mdrun[tt] produces a [TT].xvg[tt] output file that",
278         "contains projections of positions, velocities and forces onto selected",
279         "eigenvectors.[PAR]",
280         "When user-defined potential functions have been selected in the",
281         "[TT].mdp[tt] file the [TT]-table[tt] option is used to pass [TT]mdrun[tt]",
282         "a formatted table with potential functions. The file is read from",
283         "either the current directory or from the [TT]GMXLIB[tt] directory.",
284         "A number of pre-formatted tables are presented in the [TT]GMXLIB[tt] dir,",
285         "for 6-8, 6-9, 6-10, 6-11, 6-12 Lennard-Jones potentials with",
286         "normal Coulomb.",
287         "When pair interactions are present, a separate table for pair interaction",
288         "functions is read using the [TT]-tablep[tt] option.[PAR]",
289         "When tabulated bonded functions are present in the topology,",
290         "interaction functions are read using the [TT]-tableb[tt] option.",
291         "For each different tabulated interaction type the table file name is",
292         "modified in a different way: before the file extension an underscore is",
293         "appended, then a 'b' for bonds, an 'a' for angles or a 'd' for dihedrals",
294         "and finally the table number of the interaction type.[PAR]",
295         "The options [TT]-px[tt] and [TT]-pf[tt] are used for writing pull COM",
296         "coordinates and forces when pulling is selected",
297         "in the [TT].mdp[tt] file.[PAR]",
298         "With [TT]-multi[tt] or [TT]-multidir[tt], multiple systems can be ",
299         "simulated in parallel.",
300         "As many input files/directories are required as the number of systems. ",
301         "The [TT]-multidir[tt] option takes a list of directories (one for each ",
302         "system) and runs in each of them, using the input/output file names, ",
303         "such as specified by e.g. the [TT]-s[tt] option, relative to these ",
304         "directories.",
305         "With [TT]-multi[tt], the system number is appended to the run input ",
306         "and each output filename, for instance [TT]topol.tpr[tt] becomes",
307         "[TT]topol0.tpr[tt], [TT]topol1.tpr[tt] etc.",
308         "The number of nodes per system is the total number of nodes",
309         "divided by the number of systems.",
310         "One use of this option is for NMR refinement: when distance",
311         "or orientation restraints are present these can be ensemble averaged",
312         "over all the systems.[PAR]",
313         "With [TT]-replex[tt] replica exchange is attempted every given number",
314         "of steps. The number of replicas is set with the [TT]-multi[tt] or ",
315         "[TT]-multidir[tt] option, described above.",
316         "All run input files should use a different coupling temperature,",
317         "the order of the files is not important. The random seed is set with",
318         "[TT]-reseed[tt]. The velocities are scaled and neighbor searching",
319         "is performed after every exchange.[PAR]",
320         "Finally some experimental algorithms can be tested when the",
321         "appropriate options have been given. Currently under",
322         "investigation are: polarizability.",
323         "[PAR]",
324         "The option [TT]-membed[tt] does what used to be g_membed, i.e. embed",
325         "a protein into a membrane. The data file should contain the options",
326         "that where passed to g_membed before. The [TT]-mn[tt] and [TT]-mp[tt]",
327         "both apply to this as well.",
328         "[PAR]",
329         "The option [TT]-pforce[tt] is useful when you suspect a simulation",
330         "crashes due to too large forces. With this option coordinates and",
331         "forces of atoms with a force larger than a certain value will",
332         "be printed to stderr.",
333         "[PAR]",
334         "Checkpoints containing the complete state of the system are written",
335         "at regular intervals (option [TT]-cpt[tt]) to the file [TT]-cpo[tt],",
336         "unless option [TT]-cpt[tt] is set to -1.",
337         "The previous checkpoint is backed up to [TT]state_prev.cpt[tt] to",
338         "make sure that a recent state of the system is always available,",
339         "even when the simulation is terminated while writing a checkpoint.",
340         "With [TT]-cpnum[tt] all checkpoint files are kept and appended",
341         "with the step number.",
342         "A simulation can be continued by reading the full state from file",
343         "with option [TT]-cpi[tt]. This option is intelligent in the way that",
344         "if no checkpoint file is found, Gromacs just assumes a normal run and",
345         "starts from the first step of the [TT].tpr[tt] file. By default the output",
346         "will be appending to the existing output files. The checkpoint file",
347         "contains checksums of all output files, such that you will never",
348         "loose data when some output files are modified, corrupt or removed.",
349         "There are three scenarios with [TT]-cpi[tt]:[PAR]",
350         "[TT]*[tt] no files with matching names are present: new output files are written[PAR]",
351         "[TT]*[tt] all files are present with names and checksums matching those stored",
352         "in the checkpoint file: files are appended[PAR]",
353         "[TT]*[tt] otherwise no files are modified and a fatal error is generated[PAR]",
354         "With [TT]-noappend[tt] new output files are opened and the simulation",
355         "part number is added to all output file names.",
356         "Note that in all cases the checkpoint file itself is not renamed",
357         "and will be overwritten, unless its name does not match",
358         "the [TT]-cpo[tt] option.",
359         "[PAR]",
360         "With checkpointing the output is appended to previously written",
361         "output files, unless [TT]-noappend[tt] is used or none of the previous",
362         "output files are present (except for the checkpoint file).",
363         "The integrity of the files to be appended is verified using checksums",
364         "which are stored in the checkpoint file. This ensures that output can",
365         "not be mixed up or corrupted due to file appending. When only some",
366         "of the previous output files are present, a fatal error is generated",
367         "and no old output files are modified and no new output files are opened.",
368         "The result with appending will be the same as from a single run.",
369         "The contents will be binary identical, unless you use a different number",
370         "of nodes or dynamic load balancing or the FFT library uses optimizations",
371         "through timing.",
372         "[PAR]",
373         "With option [TT]-maxh[tt] a simulation is terminated and a checkpoint",
374         "file is written at the first neighbor search step where the run time",
375         "exceeds [TT]-maxh[tt]*0.99 hours.",
376         "[PAR]",
377         "When [TT]mdrun[tt] receives a TERM signal, it will set nsteps to the current",
378         "step plus one. When [TT]mdrun[tt] receives an INT signal (e.g. when ctrl+C is",
379         "pressed), it will stop after the next neighbor search step ",
380         "(with nstlist=0 at the next step).",
381         "In both cases all the usual output will be written to file.",
382         "When running with MPI, a signal to one of the [TT]mdrun[tt] processes",
383         "is sufficient, this signal should not be sent to mpirun or",
384         "the [TT]mdrun[tt] process that is the parent of the others.",
385         "[PAR]",
386         "Interactive molecular dynamics (IMD) can be activated by using at least one",
387         "of the three IMD switches: The [TT]-imdterm[tt] switch allows to terminate the",
388         "simulation from the molecular viewer (e.g. VMD). With [TT]-imdwait[tt],",
389         "[TT]mdrun[tt] pauses whenever no IMD client is connected. Pulling from the",
390         "IMD remote can be turned on by [TT]-imdpull[tt].",
391         "The port [TT]mdrun[tt] listens to can be altered by [TT]-imdport[tt].The",
392         "file pointed to by [TT]-if[tt] contains atom indices and forces if IMD",
393         "pulling is used."
394         "[PAR]",
395         "When [TT]mdrun[tt] is started with MPI, it does not run niced by default."
396     };
397     t_commrec    *cr;
398     t_filenm      fnm[] = {
399         { efTPX, NULL,      NULL,       ffREAD },
400         { efTRN, "-o",      NULL,       ffWRITE },
401         { efCOMPRESSED, "-x", NULL,     ffOPTWR },
402         { efCPT, "-cpi",    NULL,       ffOPTRD },
403         { efCPT, "-cpo",    NULL,       ffOPTWR },
404         { efSTO, "-c",      "confout",  ffWRITE },
405         { efEDR, "-e",      "ener",     ffWRITE },
406         { efLOG, "-g",      "md",       ffWRITE },
407         { efXVG, "-dhdl",   "dhdl",     ffOPTWR },
408         { efXVG, "-field",  "field",    ffOPTWR },
409         { efXVG, "-table",  "table",    ffOPTRD },
410         { efXVG, "-tabletf", "tabletf",    ffOPTRD },
411         { efXVG, "-tablep", "tablep",   ffOPTRD },
412         { efXVG, "-tableb", "table",    ffOPTRD },
413         { efTRX, "-rerun",  "rerun",    ffOPTRD },
414         { efXVG, "-tpi",    "tpi",      ffOPTWR },
415         { efXVG, "-tpid",   "tpidist",  ffOPTWR },
416         { efEDI, "-ei",     "sam",      ffOPTRD },
417         { efXVG, "-eo",     "edsam",    ffOPTWR },
418         { efXVG, "-devout", "deviatie", ffOPTWR },
419         { efXVG, "-runav",  "runaver",  ffOPTWR },
420         { efXVG, "-px",     "pullx",    ffOPTWR },
421         { efXVG, "-pf",     "pullf",    ffOPTWR },
422         { efXVG, "-ro",     "rotation", ffOPTWR },
423         { efLOG, "-ra",     "rotangles", ffOPTWR },
424         { efLOG, "-rs",     "rotslabs", ffOPTWR },
425         { efLOG, "-rt",     "rottorque", ffOPTWR },
426         { efMTX, "-mtx",    "nm",       ffOPTWR },
427         { efNDX, "-dn",     "dipole",   ffOPTWR },
428         { efRND, "-multidir", NULL,      ffOPTRDMULT},
429         { efDAT, "-membed", "membed",   ffOPTRD },
430         { efTOP, "-mp",     "membed",   ffOPTRD },
431         { efNDX, "-mn",     "membed",   ffOPTRD },
432         { efXVG, "-if",     "imdforces", ffOPTWR },
433         { efXVG, "-swap",   "swapions", ffOPTWR }
434     };
435 #define NFILE asize(fnm)
436
437     /* Command line options ! */
438     gmx_bool        bDDBondCheck  = TRUE;
439     gmx_bool        bDDBondComm   = TRUE;
440     gmx_bool        bTunePME      = TRUE;
441     gmx_bool        bTestVerlet   = FALSE;
442     gmx_bool        bVerbose      = FALSE;
443     gmx_bool        bCompact      = TRUE;
444     gmx_bool        bSepPot       = FALSE;
445     gmx_bool        bRerunVSite   = FALSE;
446     gmx_bool        bConfout      = TRUE;
447     gmx_bool        bReproducible = FALSE;
448     gmx_bool        bIMDwait      = FALSE;
449     gmx_bool        bIMDterm      = FALSE;
450     gmx_bool        bIMDpull      = FALSE;
451
452     int             npme          = -1;
453     int             nstlist       = 0;
454     int             nmultisim     = 0;
455     int             nstglobalcomm = -1;
456     int             repl_ex_nst   = 0;
457     int             repl_ex_seed  = -1;
458     int             repl_ex_nex   = 0;
459     int             nstepout      = 100;
460     int             resetstep     = -1;
461     gmx_int64_t     nsteps        = -2;   /* the value -2 means that the mdp option will be used */
462     int             imdport       = 8888; /* can be almost anything, 8888 is easy to remember */
463
464     rvec            realddxyz          = {0, 0, 0};
465     const char     *ddno_opt[ddnoNR+1] =
466     { NULL, "interleave", "pp_pme", "cartesian", NULL };
467     const char     *dddlb_opt[] =
468     { NULL, "auto", "no", "yes", NULL };
469     const char     *thread_aff_opt[threadaffNR+1] =
470     { NULL, "auto", "on", "off", NULL };
471     const char     *nbpu_opt[] =
472     { NULL, "auto", "cpu", "gpu", "gpu_cpu", NULL };
473     real            rdd                   = 0.0, rconstr = 0.0, dlb_scale = 0.8, pforce = -1;
474     char           *ddcsx                 = NULL, *ddcsy = NULL, *ddcsz = NULL;
475     real            cpt_period            = 15.0, max_hours = -1;
476     gmx_bool        bAppendFiles          = TRUE;
477     gmx_bool        bKeepAndNumCPT        = FALSE;
478     gmx_bool        bResetCountersHalfWay = FALSE;
479     output_env_t    oenv                  = NULL;
480     const char     *deviceOptions         = "";
481
482     /* Non transparent initialization of a complex gmx_hw_opt_t struct.
483      * But unfortunately we are not allowed to call a function here,
484      * since declarations follow below.
485      */
486     gmx_hw_opt_t    hw_opt = {
487         0, 0, 0, 0, threadaffSEL, 0, 0,
488         { NULL, FALSE, 0, NULL }
489     };
490
491     t_pargs         pa[] = {
492
493         { "-dd",      FALSE, etRVEC, {&realddxyz},
494           "Domain decomposition grid, 0 is optimize" },
495         { "-ddorder", FALSE, etENUM, {ddno_opt},
496           "DD node order" },
497         { "-npme",    FALSE, etINT, {&npme},
498           "Number of separate nodes to be used for PME, -1 is guess" },
499         { "-nt",      FALSE, etINT, {&hw_opt.nthreads_tot},
500           "Total number of threads to start (0 is guess)" },
501         { "-ntmpi",   FALSE, etINT, {&hw_opt.nthreads_tmpi},
502           "Number of thread-MPI threads to start (0 is guess)" },
503         { "-ntomp",   FALSE, etINT, {&hw_opt.nthreads_omp},
504           "Number of OpenMP threads per MPI process/thread to start (0 is guess)" },
505         { "-ntomp_pme", FALSE, etINT, {&hw_opt.nthreads_omp_pme},
506           "Number of OpenMP threads per MPI process/thread to start (0 is -ntomp)" },
507         { "-pin",     FALSE, etENUM, {thread_aff_opt},
508           "Fix threads (or processes) to specific cores" },
509         { "-pinoffset", FALSE, etINT, {&hw_opt.core_pinning_offset},
510           "The starting logical core number for pinning to cores; used to avoid pinning threads from different mdrun instances to the same core" },
511         { "-pinstride", FALSE, etINT, {&hw_opt.core_pinning_stride},
512           "Pinning distance in logical cores for threads, use 0 to minimize the number of threads per physical core" },
513         { "-gpu_id",  FALSE, etSTR, {&hw_opt.gpu_opt.gpu_id},
514           "List of GPU device id-s to use, specifies the per-node PP rank to GPU mapping" },
515         { "-ddcheck", FALSE, etBOOL, {&bDDBondCheck},
516           "Check for all bonded interactions with DD" },
517         { "-ddbondcomm", FALSE, etBOOL, {&bDDBondComm},
518           "HIDDENUse special bonded atom communication when [TT]-rdd[tt] > cut-off" },
519         { "-rdd",     FALSE, etREAL, {&rdd},
520           "The maximum distance for bonded interactions with DD (nm), 0 is determine from initial coordinates" },
521         { "-rcon",    FALSE, etREAL, {&rconstr},
522           "Maximum distance for P-LINCS (nm), 0 is estimate" },
523         { "-dlb",     FALSE, etENUM, {dddlb_opt},
524           "Dynamic load balancing (with DD)" },
525         { "-dds",     FALSE, etREAL, {&dlb_scale},
526           "Fraction in (0,1) by whose reciprocal the initial DD cell size will be increased in order to "
527           "provide a margin in which dynamic load balancing can act while preserving the minimum cell size." },
528         { "-ddcsx",   FALSE, etSTR, {&ddcsx},
529           "HIDDENA string containing a vector of the relative sizes in the x "
530           "direction of the corresponding DD cells. Only effective with static "
531           "load balancing." },
532         { "-ddcsy",   FALSE, etSTR, {&ddcsy},
533           "HIDDENA string containing a vector of the relative sizes in the y "
534           "direction of the corresponding DD cells. Only effective with static "
535           "load balancing." },
536         { "-ddcsz",   FALSE, etSTR, {&ddcsz},
537           "HIDDENA string containing a vector of the relative sizes in the z "
538           "direction of the corresponding DD cells. Only effective with static "
539           "load balancing." },
540         { "-gcom",    FALSE, etINT, {&nstglobalcomm},
541           "Global communication frequency" },
542         { "-nb",      FALSE, etENUM, {&nbpu_opt},
543           "Calculate non-bonded interactions on" },
544         { "-nstlist", FALSE, etINT, {&nstlist},
545           "Set nstlist when using a Verlet buffer tolerance (0 is guess)" },
546         { "-tunepme", FALSE, etBOOL, {&bTunePME},
547           "Optimize PME load between PP/PME nodes or GPU/CPU" },
548         { "-testverlet", FALSE, etBOOL, {&bTestVerlet},
549           "Test the Verlet non-bonded scheme" },
550         { "-v",       FALSE, etBOOL, {&bVerbose},
551           "Be loud and noisy" },
552         { "-compact", FALSE, etBOOL, {&bCompact},
553           "Write a compact log file" },
554         { "-seppot",  FALSE, etBOOL, {&bSepPot},
555           "Write separate V and dVdl terms for each interaction type and node to the log file(s)" },
556         { "-pforce",  FALSE, etREAL, {&pforce},
557           "Print all forces larger than this (kJ/mol nm)" },
558         { "-reprod",  FALSE, etBOOL, {&bReproducible},
559           "Try to avoid optimizations that affect binary reproducibility" },
560         { "-cpt",     FALSE, etREAL, {&cpt_period},
561           "Checkpoint interval (minutes)" },
562         { "-cpnum",   FALSE, etBOOL, {&bKeepAndNumCPT},
563           "Keep and number checkpoint files" },
564         { "-append",  FALSE, etBOOL, {&bAppendFiles},
565           "Append to previous output files when continuing from checkpoint instead of adding the simulation part number to all file names" },
566         { "-nsteps",  FALSE, etINT64, {&nsteps},
567           "Run this number of steps, overrides .mdp file option" },
568         { "-maxh",   FALSE, etREAL, {&max_hours},
569           "Terminate after 0.99 times this time (hours)" },
570         { "-multi",   FALSE, etINT, {&nmultisim},
571           "Do multiple simulations in parallel" },
572         { "-replex",  FALSE, etINT, {&repl_ex_nst},
573           "Attempt replica exchange periodically with this period (steps)" },
574         { "-nex",  FALSE, etINT, {&repl_ex_nex},
575           "Number of random exchanges to carry out each exchange interval (N^3 is one suggestion).  -nex zero or not specified gives neighbor replica exchange." },
576         { "-reseed",  FALSE, etINT, {&repl_ex_seed},
577           "Seed for replica exchange, -1 is generate a seed" },
578         { "-imdport",    FALSE, etINT, {&imdport},
579           "HIDDENIMD listening port" },
580         { "-imdwait",  FALSE, etBOOL, {&bIMDwait},
581           "HIDDENPause the simulation while no IMD client is connected" },
582         { "-imdterm",  FALSE, etBOOL, {&bIMDterm},
583           "HIDDENAllow termination of the simulation from IMD client" },
584         { "-imdpull",  FALSE, etBOOL, {&bIMDpull},
585           "HIDDENAllow pulling in the simulation from IMD client" },
586         { "-rerunvsite", FALSE, etBOOL, {&bRerunVSite},
587           "HIDDENRecalculate virtual site coordinates with [TT]-rerun[tt]" },
588         { "-confout", FALSE, etBOOL, {&bConfout},
589           "HIDDENWrite the last configuration with [TT]-c[tt] and force checkpointing at the last step" },
590         { "-stepout", FALSE, etINT, {&nstepout},
591           "HIDDENFrequency of writing the remaining wall clock time for the run" },
592         { "-resetstep", FALSE, etINT, {&resetstep},
593           "HIDDENReset cycle counters after these many time steps" },
594         { "-resethway", FALSE, etBOOL, {&bResetCountersHalfWay},
595           "HIDDENReset the cycle counters after half the number of steps or halfway [TT]-maxh[tt]" }
596     };
597     unsigned long   Flags, PCA_Flags;
598     ivec            ddxyz;
599     int             dd_node_order;
600     gmx_bool        bAddPart;
601     FILE           *fplog, *fpmulti;
602     int             sim_part, sim_part_fn;
603     const char     *part_suffix = ".part";
604     char            suffix[STRLEN];
605     int             rc;
606     char          **multidir = NULL;
607
608
609     cr = init_commrec();
610
611     PCA_Flags = (PCA_CAN_SET_DEFFNM | (MASTER(cr) ? 0 : PCA_QUIET));
612
613     /* Comment this in to do fexist calls only on master
614      * works not with rerun or tables at the moment
615      * also comment out the version of init_forcerec in md.c
616      * with NULL instead of opt2fn
617      */
618     /*
619        if (!MASTER(cr))
620        {
621        PCA_Flags |= PCA_NOT_READ_NODE;
622        }
623      */
624
625     if (!parse_common_args(&argc, argv, PCA_Flags, NFILE, fnm, asize(pa), pa,
626                            asize(desc), desc, 0, NULL, &oenv))
627     {
628         return 0;
629     }
630
631
632     /* we set these early because they might be used in init_multisystem()
633        Note that there is the potential for npme>nnodes until the number of
634        threads is set later on, if there's thread parallelization. That shouldn't
635        lead to problems. */
636     dd_node_order = nenum(ddno_opt);
637     cr->npmenodes = npme;
638
639     hw_opt.thread_affinity = nenum(thread_aff_opt);
640
641     /* now check the -multi and -multidir option */
642     if (opt2bSet("-multidir", NFILE, fnm))
643     {
644         if (nmultisim > 0)
645         {
646             gmx_fatal(FARGS, "mdrun -multi and -multidir options are mutually exclusive.");
647         }
648         nmultisim = opt2fns(&multidir, "-multidir", NFILE, fnm);
649     }
650
651
652     if (repl_ex_nst != 0 && nmultisim < 2)
653     {
654         gmx_fatal(FARGS, "Need at least two replicas for replica exchange (option -multi)");
655     }
656
657     if (repl_ex_nex < 0)
658     {
659         gmx_fatal(FARGS, "Replica exchange number of exchanges needs to be positive");
660     }
661
662     if (nmultisim > 1)
663     {
664 #ifndef GMX_THREAD_MPI
665         gmx_bool bParFn = (multidir == NULL);
666         init_multisystem(cr, nmultisim, multidir, NFILE, fnm, bParFn);
667 #else
668         gmx_fatal(FARGS, "mdrun -multi is not supported with the thread library. "
669                   "Please compile GROMACS with MPI support");
670 #endif
671     }
672
673     bAddPart = !bAppendFiles;
674
675     /* Check if there is ANY checkpoint file available */
676     sim_part    = 1;
677     sim_part_fn = sim_part;
678     if (opt2bSet("-cpi", NFILE, fnm))
679     {
680         if (bSepPot && bAppendFiles)
681         {
682             gmx_fatal(FARGS, "Output file appending is not supported with -seppot");
683         }
684
685         bAppendFiles =
686             read_checkpoint_simulation_part(opt2fn_master("-cpi", NFILE,
687                                                           fnm, cr),
688                                             &sim_part_fn, NULL, cr,
689                                             bAppendFiles, NFILE, fnm,
690                                             part_suffix, &bAddPart);
691         if (sim_part_fn == 0 && MULTIMASTER(cr))
692         {
693             fprintf(stdout, "No previous checkpoint file present, assuming this is a new run.\n");
694         }
695         else
696         {
697             sim_part = sim_part_fn + 1;
698         }
699
700         if (MULTISIM(cr) && MASTER(cr))
701         {
702             if (MULTIMASTER(cr))
703             {
704                 /* Log file is not yet available, so if there's a
705                  * problem we can only write to stderr. */
706                 fpmulti = stderr;
707             }
708             else
709             {
710                 fpmulti = NULL;
711             }
712             check_multi_int(fpmulti, cr->ms, sim_part, "simulation part", TRUE);
713         }
714     }
715     else
716     {
717         bAppendFiles = FALSE;
718     }
719
720     if (!bAppendFiles)
721     {
722         sim_part_fn = sim_part;
723     }
724
725     if (bAddPart)
726     {
727         /* Rename all output files (except checkpoint files) */
728         /* create new part name first (zero-filled) */
729         sprintf(suffix, "%s%04d", part_suffix, sim_part_fn);
730
731         add_suffix_to_output_names(fnm, NFILE, suffix);
732         if (MULTIMASTER(cr))
733         {
734             fprintf(stdout, "Checkpoint file is from part %d, new output files will be suffixed '%s'.\n", sim_part-1, suffix);
735         }
736     }
737
738     Flags = opt2bSet("-rerun", NFILE, fnm) ? MD_RERUN : 0;
739     Flags = Flags | (bSepPot       ? MD_SEPPOT       : 0);
740     Flags = Flags | (bDDBondCheck  ? MD_DDBONDCHECK  : 0);
741     Flags = Flags | (bDDBondComm   ? MD_DDBONDCOMM   : 0);
742     Flags = Flags | (bTunePME      ? MD_TUNEPME      : 0);
743     Flags = Flags | (bTestVerlet   ? MD_TESTVERLET   : 0);
744     Flags = Flags | (bConfout      ? MD_CONFOUT      : 0);
745     Flags = Flags | (bRerunVSite   ? MD_RERUN_VSITE  : 0);
746     Flags = Flags | (bReproducible ? MD_REPRODUCIBLE : 0);
747     Flags = Flags | (bAppendFiles  ? MD_APPENDFILES  : 0);
748     Flags = Flags | (opt2parg_bSet("-append", asize(pa), pa) ? MD_APPENDFILESSET : 0);
749     Flags = Flags | (bKeepAndNumCPT ? MD_KEEPANDNUMCPT : 0);
750     Flags = Flags | (sim_part > 1    ? MD_STARTFROMCPT : 0);
751     Flags = Flags | (bResetCountersHalfWay ? MD_RESETCOUNTERSHALFWAY : 0);
752     Flags = Flags | (bIMDwait      ? MD_IMDWAIT      : 0);
753     Flags = Flags | (bIMDterm      ? MD_IMDTERM      : 0);
754     Flags = Flags | (bIMDpull      ? MD_IMDPULL      : 0);
755
756     /* We postpone opening the log file if we are appending, so we can
757        first truncate the old log file and append to the correct position
758        there instead.  */
759     if ((MASTER(cr) || bSepPot) && !bAppendFiles)
760     {
761         gmx_log_open(ftp2fn(efLOG, NFILE, fnm), cr,
762                      !bSepPot, Flags & MD_APPENDFILES, &fplog);
763         please_cite(fplog, "Hess2008b");
764         please_cite(fplog, "Spoel2005a");
765         please_cite(fplog, "Lindahl2001a");
766         please_cite(fplog, "Berendsen95a");
767     }
768     else if (!MASTER(cr) && bSepPot)
769     {
770         gmx_log_open(ftp2fn(efLOG, NFILE, fnm), cr, !bSepPot, Flags, &fplog);
771     }
772     else
773     {
774         fplog = NULL;
775     }
776
777     ddxyz[XX] = (int)(realddxyz[XX] + 0.5);
778     ddxyz[YY] = (int)(realddxyz[YY] + 0.5);
779     ddxyz[ZZ] = (int)(realddxyz[ZZ] + 0.5);
780
781     rc = mdrunner(&hw_opt, fplog, cr, NFILE, fnm, oenv, bVerbose, bCompact,
782                   nstglobalcomm, ddxyz, dd_node_order, rdd, rconstr,
783                   dddlb_opt[0], dlb_scale, ddcsx, ddcsy, ddcsz,
784                   nbpu_opt[0], nstlist,
785                   nsteps, nstepout, resetstep,
786                   nmultisim, repl_ex_nst, repl_ex_nex, repl_ex_seed,
787                   pforce, cpt_period, max_hours, deviceOptions, imdport, Flags);
788
789     /* Log file has to be closed in mdrunner if we are appending to it
790        (fplog not set here) */
791     if (MASTER(cr) && !bAppendFiles)
792     {
793         gmx_log_close(fplog);
794     }
795
796     return rc;
797 }