1d011082be2a359a5c38c0576bea30aeccfcd35a
[alexxy/gromacs.git] / src / kernel / mdrun.c
1 /*  -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  * 
4  *                This source code is part of
5  * 
6  *                 G   R   O   M   A   C   S
7  * 
8  *          GROningen MAchine for Chemical Simulations
9  * 
10  *                        VERSION 3.2.0
11  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13  * Copyright (c) 2001-2004, The GROMACS development team,
14  * check out http://www.gromacs.org for more information.
15
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  * 
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  * 
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  * 
31  * For more info, check our website at http://www.gromacs.org
32  * 
33  * And Hey:
34  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include "typedefs.h"
41 #include "macros.h"
42 #include "copyrite.h"
43 #include "main.h"
44 #include "statutil.h"
45 #include "smalloc.h"
46 #include "futil.h"
47 #include "smalloc.h"
48 #include "edsam.h"
49 #include "mdrun.h"
50 #include "xmdrun.h"
51 #include "checkpoint.h"
52 #ifdef GMX_THREAD_MPI
53 #include "thread_mpi.h"
54 #endif
55
56 /* afm stuf */
57 #include "pull.h"
58
59 int main(int argc,char *argv[])
60 {
61   const char *desc[] = {
62  #ifdef GMX_OPENMM
63     "This is an experimental release of GROMACS for accelerated",
64         "Molecular Dynamics simulations on GPU processors. Support is provided",
65         "by the OpenMM library (https://simtk.org/home/openmm).[PAR]",
66         "*Warning*[BR]",
67         "This release is targeted at developers and advanced users and",
68         "care should be taken before production use. The following should be",
69         "noted before using the program:[PAR]",
70         " * The current release runs only on modern nVidia GPU hardware with CUDA support.",
71         "Make sure that the necessary CUDA drivers and libraries for your operating system",
72         "are already installed. The CUDA SDK also should be installed in order to compile",
73         "the program from source (http://www.nvidia.com/object/cuda_home.html).[PAR]",
74         " * Multiple GPU cards are not supported.[PAR]",
75         " * Only a small subset of the GROMACS features and options are supported on the GPUs.",
76         "See below for a detailed list.[PAR]",
77         " * Consumer level GPU cards are known to often have problems with faulty memory.",
78         "It is recommended that a full memory check of the cards is done at least once",
79         "(for example, using the memtest=full option).",
80         "A partial memory check (for example, memtest=15) before and",
81         "after the simulation run would help spot",
82         "problems resulting from processor overheating.[PAR]",
83         " * The maximum size of the simulated systems depends on the available",
84         "GPU memory,for example, a GTX280 with 1GB memory has been tested with systems",
85         "of up to about 100,000 atoms.[PAR]",
86         " * In order to take a full advantage of the GPU platform features, many algorithms",
87         "have been implemented in a very different way than they are on the CPUs.",
88         "Therefore numercal correspondence between properties of the state of",
89         "simulated systems should not be expected. Moreover, the values will likely vary",
90         "when simulations are done on different GPU hardware.[PAR]",
91         " * Frequent retrieval of system state information such as",
92         "trajectory coordinates and energies can greatly influence the performance",
93         "of the program due to slow CPU<->GPU memory transfer speed.[PAR]",
94         " * MD algorithms are complex, and although the Gromacs code is highly tuned for them,",
95         "they often do not translate very well onto the streaming architetures.",
96         "Realistic expectations about the achievable speed-up from test with GTX280:",
97         "For small protein systems in implicit solvent using all-vs-all kernels the acceleration",
98         "can be as high as 20 times, but in most other setups involving cutoffs and PME the",
99         "acceleration is usually only ~4 times relative to a 3GHz CPU.[PAR]",
100         "Supported features:[PAR]",
101         " * Integrators: md/md-vv/md-vv-avek, sd/sd1 and bd.\n",
102         " * Long-range interactions (option coulombtype): Reaction-Field, Ewald, PME, and cut-off (for Implicit Solvent only)\n",
103         " * Temperature control: Supported only with the md/md-vv/md-vv-avek, sd/sd1 and bd integrators.\n",
104         " * Pressure control: Supported.\n",
105         " * Implicit solvent: Supported.\n",
106         "A detailed description can be found on the GROMACS website:\n",
107         "http://www.gromacs.org/gpu[PAR]",
108 /* From the original mdrun documentaion */
109     "The [TT]mdrun[tt] program reads the run input file ([TT]-s[tt])",
110     "and distributes the topology over nodes if needed.",
111     "[TT]mdrun[tt] produces at least four output files.",
112     "A single log file ([TT]-g[tt]) is written, unless the option",
113     "[TT]-seppot[tt] is used, in which case each node writes a log file.",
114     "The trajectory file ([TT]-o[tt]), contains coordinates, velocities and",
115     "optionally forces.",
116     "The structure file ([TT]-c[tt]) contains the coordinates and",
117     "velocities of the last step.",
118     "The energy file ([TT]-e[tt]) contains energies, the temperature,",
119     "pressure, etc, a lot of these things are also printed in the log file.",
120     "Optionally coordinates can be written to a compressed trajectory file",
121     "([TT]-x[tt]).[PAR]",
122 /* openmm specific information */
123         "Usage with OpenMM:[BR]",
124         "[TT]mdrun -device \"OpenMM:platform=Cuda,memtest=15,deviceid=0,force-device=no\"[tt][PAR]",
125         "Options:[PAR]",
126         "      [TT]platform[tt] = Cuda\t\t:\tThe only available value. OpenCL support will be available in future.\n",
127         "      [TT]memtest[tt] = 15\t\t:\tRun a partial, random GPU memory test for the given amount of seconds. A full test",
128         "(recommended!) can be run with \"memtest=full\". Memory testing can be disabled with \"memtest=off\".\n",
129         "      [TT]deviceid[tt] = 0\t\t:\tSpecify the target device when multiple cards are present.",
130         "Only one card can be used at any given time though.\n",
131         "      [TT]force-device[tt] = no\t\t:\tIf set to \"yes\" [TT]mdrun[tt]  will be forced to execute on",
132         "hardware that is not officially supported. GPU acceleration can also be achieved on older",
133         "but Cuda capable cards, although the simulation might be too slow, and the memory limits too strict.",
134 #else
135     "The [TT]mdrun[tt] program is the main computational chemistry engine",
136     "within GROMACS. Obviously, it performs Molecular Dynamics simulations,",
137     "but it can also perform Stochastic Dynamics, Energy Minimization,",
138     "test particle insertion or (re)calculation of energies.",
139     "Normal mode analysis is another option. In this case [TT]mdrun[tt]",
140     "builds a Hessian matrix from single conformation.",
141     "For usual Normal Modes-like calculations, make sure that",
142     "the structure provided is properly energy-minimized.",
143     "The generated matrix can be diagonalized by [TT]g_nmeig[tt].[PAR]",
144     "The [TT]mdrun[tt] program reads the run input file ([TT]-s[tt])",
145     "and distributes the topology over nodes if needed.",
146     "[TT]mdrun[tt] produces at least four output files.",
147     "A single log file ([TT]-g[tt]) is written, unless the option",
148     "[TT]-seppot[tt] is used, in which case each node writes a log file.",
149     "The trajectory file ([TT]-o[tt]), contains coordinates, velocities and",
150     "optionally forces.",
151     "The structure file ([TT]-c[tt]) contains the coordinates and",
152     "velocities of the last step.",
153     "The energy file ([TT]-e[tt]) contains energies, the temperature,",
154     "pressure, etc, a lot of these things are also printed in the log file.",
155     "Optionally coordinates can be written to a compressed trajectory file",
156     "([TT]-x[tt]).[PAR]",
157     "The option [TT]-dhdl[tt] is only used when free energy calculation is",
158     "turned on.[PAR]",
159     "When [TT]mdrun[tt] is started using MPI with more than 1 node, parallelization",
160     "is used. By default domain decomposition is used, unless the [TT]-pd[tt]",
161     "option is set, which selects particle decomposition.[PAR]",
162     "With domain decomposition, the spatial decomposition can be set",
163     "with option [TT]-dd[tt]. By default [TT]mdrun[tt] selects a good decomposition.",
164     "The user only needs to change this when the system is very inhomogeneous.",
165     "Dynamic load balancing is set with the option [TT]-dlb[tt],",
166     "which can give a significant performance improvement,",
167     "especially for inhomogeneous systems. The only disadvantage of",
168     "dynamic load balancing is that runs are no longer binary reproducible,",
169     "but in most cases this is not important.",
170     "By default the dynamic load balancing is automatically turned on",
171     "when the measured performance loss due to load imbalance is 5% or more.",
172     "At low parallelization these are the only important options",
173     "for domain decomposition.",
174     "At high parallelization the options in the next two sections",
175     "could be important for increasing the performace.",
176     "[PAR]",
177     "When PME is used with domain decomposition, separate nodes can",
178     "be assigned to do only the PME mesh calculation;",
179     "this is computationally more efficient starting at about 12 nodes.",
180     "The number of PME nodes is set with option [TT]-npme[tt],",
181     "this can not be more than half of the nodes.",
182     "By default [TT]mdrun[tt] makes a guess for the number of PME",
183     "nodes when the number of nodes is larger than 11 or performance wise",
184     "not compatible with the PME grid x dimension.",
185     "But the user should optimize npme. Performance statistics on this issue",
186     "are written at the end of the log file.",
187     "For good load balancing at high parallelization, the PME grid x and y",
188     "dimensions should be divisible by the number of PME nodes",
189     "(the simulation will run correctly also when this is not the case).",
190     "[PAR]",
191     "This section lists all options that affect the domain decomposition.",
192     "[PAR]",
193     "Option [TT]-rdd[tt] can be used to set the required maximum distance",
194     "for inter charge-group bonded interactions.",
195     "Communication for two-body bonded interactions below the non-bonded",
196     "cut-off distance always comes for free with the non-bonded communication.",
197     "Atoms beyond the non-bonded cut-off are only communicated when they have",
198     "missing bonded interactions; this means that the extra cost is minor",
199     "and nearly indepedent of the value of [TT]-rdd[tt].",
200     "With dynamic load balancing option [TT]-rdd[tt] also sets",
201     "the lower limit for the domain decomposition cell sizes.",
202     "By default [TT]-rdd[tt] is determined by [TT]mdrun[tt] based on",
203     "the initial coordinates. The chosen value will be a balance",
204     "between interaction range and communication cost.",
205     "[PAR]",
206     "When inter charge-group bonded interactions are beyond",
207     "the bonded cut-off distance, [TT]mdrun[tt] terminates with an error message.",
208     "For pair interactions and tabulated bonds",
209     "that do not generate exclusions, this check can be turned off",
210     "with the option [TT]-noddcheck[tt].",
211     "[PAR]",
212     "When constraints are present, option [TT]-rcon[tt] influences",
213     "the cell size limit as well.",
214     "Atoms connected by NC constraints, where NC is the LINCS order plus 1,",
215     "should not be beyond the smallest cell size. A error message is",
216     "generated when this happens and the user should change the decomposition",
217     "or decrease the LINCS order and increase the number of LINCS iterations.",
218     "By default [TT]mdrun[tt] estimates the minimum cell size required for P-LINCS",
219     "in a conservative fashion. For high parallelization it can be useful",
220     "to set the distance required for P-LINCS with the option [TT]-rcon[tt].",
221     "[PAR]",
222     "The [TT]-dds[tt] option sets the minimum allowed x, y and/or z scaling",
223     "of the cells with dynamic load balancing. [TT]mdrun[tt] will ensure that",
224     "the cells can scale down by at least this factor. This option is used",
225     "for the automated spatial decomposition (when not using [TT]-dd[tt])",
226     "as well as for determining the number of grid pulses, which in turn",
227     "sets the minimum allowed cell size. Under certain circumstances",
228     "the value of [TT]-dds[tt] might need to be adjusted to account for",
229     "high or low spatial inhomogeneity of the system.",
230     "[PAR]",
231     "The option [TT]-gcom[tt] can be used to only do global communication",
232     "every n steps.",
233     "This can improve performance for highly parallel simulations",
234     "where this global communication step becomes the bottleneck.",
235     "For a global thermostat and/or barostat the temperature",
236     "and/or pressure will also only be updated every [TT]-gcom[tt] steps.",
237     "By default it is set to the minimum of nstcalcenergy and nstlist.[PAR]",
238     "With [TT]-rerun[tt] an input trajectory can be given for which ",
239     "forces and energies will be (re)calculated. Neighbor searching will be",
240     "performed for every frame, unless [TT]nstlist[tt] is zero",
241     "(see the [TT].mdp[tt] file).[PAR]",
242     "ED (essential dynamics) sampling is switched on by using the [TT]-ei[tt]",
243     "flag followed by an [TT].edi[tt] file.",
244     "The [TT].edi[tt] file can be produced using options in the essdyn",
245     "menu of the WHAT IF program. [TT]mdrun[tt] produces a [TT].edo[tt] file that",
246     "contains projections of positions, velocities and forces onto selected",
247     "eigenvectors.[PAR]",
248     "When user-defined potential functions have been selected in the",
249     "[TT].mdp[tt] file the [TT]-table[tt] option is used to pass [TT]mdrun[tt]",
250     "a formatted table with potential functions. The file is read from",
251     "either the current directory or from the [TT]GMXLIB[tt] directory.",
252     "A number of pre-formatted tables are presented in the [TT]GMXLIB[tt] dir,",
253     "for 6-8, 6-9, 6-10, 6-11, 6-12 Lennard-Jones potentials with",
254     "normal Coulomb.",
255     "When pair interactions are present, a separate table for pair interaction",
256     "functions is read using the [TT]-tablep[tt] option.[PAR]",
257     "When tabulated bonded functions are present in the topology,",
258     "interaction functions are read using the [TT]-tableb[tt] option.",
259     "For each different tabulated interaction type the table file name is",
260     "modified in a different way: before the file extension an underscore is",
261     "appended, then a 'b' for bonds, an 'a' for angles or a 'd' for dihedrals",
262     "and finally the table number of the interaction type.[PAR]",
263     "The options [TT]-px[tt] and [TT]-pf[tt] are used for writing pull COM",
264     "coordinates and forces when pulling is selected",
265     "in the [TT].mdp[tt] file.[PAR]",
266     "With [TT]-multi[tt] or [TT]-multidir[tt], multiple systems can be ",
267     "simulated in parallel.",
268     "As many input files/directories are required as the number of systems. ",
269     "The [TT]-multidir[tt] option takes a list of directories (one for each ",
270     "system) and runs in each of them, using the input/output file names, ",
271     "such as specified by e.g. the [TT]-s[tt] option, relative to these ",
272     "directories.",
273     "With [TT]-multi[tt], the system number is appended to the run input ",
274     "and each output filename, for instance [TT]topol.tpr[tt] becomes",
275     "[TT]topol0.tpr[tt], [TT]topol1.tpr[tt] etc.",
276     "The number of nodes per system is the total number of nodes",
277     "divided by the number of systems.",
278     "One use of this option is for NMR refinement: when distance",
279     "or orientation restraints are present these can be ensemble averaged",
280     "over all the systems.[PAR]",
281     "With [TT]-replex[tt] replica exchange is attempted every given number",
282     "of steps. The number of replicas is set with the [TT]-multi[tt] or ",
283     "[TT]-multidir[tt] option, described above.",
284     "All run input files should use a different coupling temperature,",
285     "the order of the files is not important. The random seed is set with",
286     "[TT]-reseed[tt]. The velocities are scaled and neighbor searching",
287     "is performed after every exchange.[PAR]",
288     "Finally some experimental algorithms can be tested when the",
289     "appropriate options have been given. Currently under",
290     "investigation are: polarizability and X-ray bombardments.",
291     "[PAR]",
292     "The option [TT]-membed[tt] does what used to be g_membed, i.e. embed",
293     "a protein into a membrane. The data file should contain the options",
294     "that where passed to g_membed before. The [TT]-mn[tt] and [TT]-mp[tt]",
295     "both apply to this as well.",
296     "[PAR]",
297     "The option [TT]-pforce[tt] is useful when you suspect a simulation",
298     "crashes due to too large forces. With this option coordinates and",
299     "forces of atoms with a force larger than a certain value will",
300     "be printed to stderr.",
301     "[PAR]",
302     "Checkpoints containing the complete state of the system are written",
303     "at regular intervals (option [TT]-cpt[tt]) to the file [TT]-cpo[tt],",
304     "unless option [TT]-cpt[tt] is set to -1.",
305     "The previous checkpoint is backed up to [TT]state_prev.cpt[tt] to",
306     "make sure that a recent state of the system is always available,",
307     "even when the simulation is terminated while writing a checkpoint.",
308     "With [TT]-cpnum[tt] all checkpoint files are kept and appended",
309     "with the step number.",
310     "A simulation can be continued by reading the full state from file",
311     "with option [TT]-cpi[tt]. This option is intelligent in the way that",
312     "if no checkpoint file is found, Gromacs just assumes a normal run and",
313     "starts from the first step of the [TT].tpr[tt] file. By default the output",
314     "will be appending to the existing output files. The checkpoint file",
315     "contains checksums of all output files, such that you will never",
316     "loose data when some output files are modified, corrupt or removed.",
317     "There are three scenarios with [TT]-cpi[tt]:[PAR]",
318     "[TT]*[tt] no files with matching names are present: new output files are written[PAR]",
319     "[TT]*[tt] all files are present with names and checksums matching those stored",
320     "in the checkpoint file: files are appended[PAR]",
321     "[TT]*[tt] otherwise no files are modified and a fatal error is generated[PAR]",
322     "With [TT]-noappend[tt] new output files are opened and the simulation",
323     "part number is added to all output file names.",
324     "Note that in all cases the checkpoint file itself is not renamed",
325     "and will be overwritten, unless its name does not match",
326     "the [TT]-cpo[tt] option.",
327     "[PAR]",
328     "With checkpointing the output is appended to previously written",
329     "output files, unless [TT]-noappend[tt] is used or none of the previous",
330     "output files are present (except for the checkpoint file).",
331     "The integrity of the files to be appended is verified using checksums",
332     "which are stored in the checkpoint file. This ensures that output can",
333     "not be mixed up or corrupted due to file appending. When only some",
334     "of the previous output files are present, a fatal error is generated",
335     "and no old output files are modified and no new output files are opened.",
336     "The result with appending will be the same as from a single run.",
337     "The contents will be binary identical, unless you use a different number",
338     "of nodes or dynamic load balancing or the FFT library uses optimizations",
339     "through timing.",
340     "[PAR]",
341     "With option [TT]-maxh[tt] a simulation is terminated and a checkpoint",
342     "file is written at the first neighbor search step where the run time",
343     "exceeds [TT]-maxh[tt]*0.99 hours.",
344     "[PAR]",
345     "When [TT]mdrun[tt] receives a TERM signal, it will set nsteps to the current",
346     "step plus one. When [TT]mdrun[tt] receives an INT signal (e.g. when ctrl+C is",
347     "pressed), it will stop after the next neighbor search step ",
348     "(with nstlist=0 at the next step).",
349     "In both cases all the usual output will be written to file.",
350     "When running with MPI, a signal to one of the [TT]mdrun[tt] processes",
351     "is sufficient, this signal should not be sent to mpirun or",
352     "the [TT]mdrun[tt] process that is the parent of the others.",
353     "[PAR]",
354     "When [TT]mdrun[tt] is started with MPI, it does not run niced by default."
355 #endif
356   };
357   t_commrec    *cr;
358   t_filenm fnm[] = {
359     { efTPX, NULL,      NULL,       ffREAD },
360     { efTRN, "-o",      NULL,       ffWRITE },
361     { efXTC, "-x",      NULL,       ffOPTWR },
362     { efCPT, "-cpi",    NULL,       ffOPTRD },
363     { efCPT, "-cpo",    NULL,       ffOPTWR },
364     { efSTO, "-c",      "confout",  ffWRITE },
365     { efEDR, "-e",      "ener",     ffWRITE },
366     { efLOG, "-g",      "md",       ffWRITE },
367     { efXVG, "-dhdl",   "dhdl",     ffOPTWR },
368     { efXVG, "-field",  "field",    ffOPTWR },
369     { efXVG, "-table",  "table",    ffOPTRD },
370     { efXVG, "-tabletf", "tabletf",    ffOPTRD },
371     { efXVG, "-tablep", "tablep",   ffOPTRD },
372     { efXVG, "-tableb", "table",    ffOPTRD },
373     { efTRX, "-rerun",  "rerun",    ffOPTRD },
374     { efXVG, "-tpi",    "tpi",      ffOPTWR },
375     { efXVG, "-tpid",   "tpidist",  ffOPTWR },
376     { efEDI, "-ei",     "sam",      ffOPTRD },
377     { efEDO, "-eo",     "sam",      ffOPTWR },
378     { efGCT, "-j",      "wham",     ffOPTRD },
379     { efGCT, "-jo",     "bam",      ffOPTWR },
380     { efXVG, "-ffout",  "gct",      ffOPTWR },
381     { efXVG, "-devout", "deviatie", ffOPTWR },
382     { efXVG, "-runav",  "runaver",  ffOPTWR },
383     { efXVG, "-px",     "pullx",    ffOPTWR },
384     { efXVG, "-pf",     "pullf",    ffOPTWR },
385     { efXVG, "-ro",     "rotation", ffOPTWR },
386     { efLOG, "-ra",     "rotangles",ffOPTWR },
387     { efLOG, "-rs",     "rotslabs", ffOPTWR },
388     { efLOG, "-rt",     "rottorque",ffOPTWR },
389     { efMTX, "-mtx",    "nm",       ffOPTWR },
390     { efNDX, "-dn",     "dipole",   ffOPTWR },
391     { efRND, "-multidir",NULL,      ffOPTRDMULT},
392     { efDAT, "-membed", "membed",   ffOPTRD },
393     { efTOP, "-mp",     "membed",   ffOPTRD },
394     { efNDX, "-mn",     "membed",   ffOPTRD }
395   };
396 #define NFILE asize(fnm)
397
398   /* Command line options ! */
399   gmx_bool bCart        = FALSE;
400   gmx_bool bPPPME       = FALSE;
401   gmx_bool bPartDec     = FALSE;
402   gmx_bool bDDBondCheck = TRUE;
403   gmx_bool bDDBondComm  = TRUE;
404   gmx_bool bVerbose     = FALSE;
405   gmx_bool bCompact     = TRUE;
406   gmx_bool bSepPot      = FALSE;
407   gmx_bool bRerunVSite  = FALSE;
408   gmx_bool bIonize      = FALSE;
409   gmx_bool bConfout     = TRUE;
410   gmx_bool bReproducible = FALSE;
411     
412   int  npme=-1;
413   int  nmultisim=0;
414   int  nstglobalcomm=-1;
415   int  repl_ex_nst=0;
416   int  repl_ex_seed=-1;
417   int  repl_ex_nex=0;
418   int  nstepout=100;
419   int  nthreads=0; /* set to determine # of threads automatically */
420   int  resetstep=-1;
421   
422   rvec realddxyz={0,0,0};
423   const char *ddno_opt[ddnoNR+1] =
424     { NULL, "interleave", "pp_pme", "cartesian", NULL };
425     const char *dddlb_opt[] =
426     { NULL, "auto", "no", "yes", NULL };
427   real rdd=0.0,rconstr=0.0,dlb_scale=0.8,pforce=-1;
428   char *ddcsx=NULL,*ddcsy=NULL,*ddcsz=NULL;
429   real cpt_period=15.0,max_hours=-1;
430   gmx_bool bAppendFiles=TRUE;
431   gmx_bool bKeepAndNumCPT=FALSE;
432   gmx_bool bResetCountersHalfWay=FALSE;
433   output_env_t oenv=NULL;
434   const char *deviceOptions = "";
435
436   t_pargs pa[] = {
437
438     { "-pd",      FALSE, etBOOL,{&bPartDec},
439       "Use particle decompostion" },
440     { "-dd",      FALSE, etRVEC,{&realddxyz},
441       "Domain decomposition grid, 0 is optimize" },
442 #ifdef GMX_THREAD_MPI
443     { "-nt",      FALSE, etINT, {&nthreads},
444       "Number of threads to start (0 is guess)" },
445 #endif
446     { "-npme",    FALSE, etINT, {&npme},
447       "Number of separate nodes to be used for PME, -1 is guess" },
448     { "-ddorder", FALSE, etENUM, {ddno_opt},
449       "DD node order" },
450     { "-ddcheck", FALSE, etBOOL, {&bDDBondCheck},
451       "Check for all bonded interactions with DD" },
452     { "-ddbondcomm", FALSE, etBOOL, {&bDDBondComm},
453       "HIDDENUse special bonded atom communication when [TT]-rdd[tt] > cut-off" },
454     { "-rdd",     FALSE, etREAL, {&rdd},
455       "The maximum distance for bonded interactions with DD (nm), 0 is determine from initial coordinates" },
456     { "-rcon",    FALSE, etREAL, {&rconstr},
457       "Maximum distance for P-LINCS (nm), 0 is estimate" },
458     { "-dlb",     FALSE, etENUM, {dddlb_opt},
459       "Dynamic load balancing (with DD)" },
460     { "-dds",     FALSE, etREAL, {&dlb_scale},
461       "Minimum allowed dlb scaling of the DD cell size" },
462     { "-ddcsx",   FALSE, etSTR, {&ddcsx},
463       "HIDDENThe DD cell sizes in x" },
464     { "-ddcsy",   FALSE, etSTR, {&ddcsy},
465       "HIDDENThe DD cell sizes in y" },
466     { "-ddcsz",   FALSE, etSTR, {&ddcsz},
467       "HIDDENThe DD cell sizes in z" },
468     { "-gcom",    FALSE, etINT,{&nstglobalcomm},
469       "Global communication frequency" },
470     { "-v",       FALSE, etBOOL,{&bVerbose},  
471       "Be loud and noisy" },
472     { "-compact", FALSE, etBOOL,{&bCompact},  
473       "Write a compact log file" },
474     { "-seppot",  FALSE, etBOOL, {&bSepPot},
475       "Write separate V and dVdl terms for each interaction type and node to the log file(s)" },
476     { "-pforce",  FALSE, etREAL, {&pforce},
477       "Print all forces larger than this (kJ/mol nm)" },
478     { "-reprod",  FALSE, etBOOL,{&bReproducible},  
479       "Try to avoid optimizations that affect binary reproducibility" },
480     { "-cpt",     FALSE, etREAL, {&cpt_period},
481       "Checkpoint interval (minutes)" },
482     { "-cpnum",   FALSE, etBOOL, {&bKeepAndNumCPT},
483       "Keep and number checkpoint files" },
484     { "-append",  FALSE, etBOOL, {&bAppendFiles},
485       "Append to previous output files when continuing from checkpoint instead of adding the simulation part number to all file names" },
486     { "-maxh",   FALSE, etREAL, {&max_hours},
487       "Terminate after 0.99 times this time (hours)" },
488     { "-multi",   FALSE, etINT,{&nmultisim}, 
489       "Do multiple simulations in parallel" },
490     { "-replex",  FALSE, etINT, {&repl_ex_nst}, 
491       "Attempt replica exchange periodically with this period (steps)" },
492     { "-nex",  FALSE, etINT, {&repl_ex_nex},
493       "Number of random exchanges to carry out each exchange interval (N^3 is one suggestion).  -nex zero or not specified gives neighbor replica exchange." },
494     { "-reseed",  FALSE, etINT, {&repl_ex_seed}, 
495       "Seed for replica exchange, -1 is generate a seed" },
496     { "-rerunvsite", FALSE, etBOOL, {&bRerunVSite},
497       "HIDDENRecalculate virtual site coordinates with [TT]-rerun[tt]" },
498     { "-ionize",  FALSE, etBOOL,{&bIonize},
499       "Do a simulation including the effect of an X-Ray bombardment on your system" },
500     { "-confout", FALSE, etBOOL, {&bConfout},
501       "HIDDENWrite the last configuration with [TT]-c[tt] and force checkpointing at the last step" },
502     { "-stepout", FALSE, etINT, {&nstepout},
503       "HIDDENFrequency of writing the remaining runtime" },
504     { "-resetstep", FALSE, etINT, {&resetstep},
505       "HIDDENReset cycle counters after these many time steps" },
506     { "-resethway", FALSE, etBOOL, {&bResetCountersHalfWay},
507       "HIDDENReset the cycle counters after half the number of steps or halfway [TT]-maxh[tt]" }
508 #ifdef GMX_OPENMM
509     ,
510     { "-device",  FALSE, etSTR, {&deviceOptions},
511       "Device option string" }
512 #endif
513   };
514   gmx_edsam_t  ed;
515   unsigned long Flags, PCA_Flags;
516   ivec     ddxyz;
517   int      dd_node_order;
518   gmx_bool     bAddPart;
519   FILE     *fplog,*fptest;
520   int      sim_part,sim_part_fn;
521   const char *part_suffix=".part";
522   char     suffix[STRLEN];
523   int      rc;
524   char **multidir=NULL;
525
526
527   cr = init_par(&argc,&argv);
528
529   if (MASTER(cr))
530     CopyRight(stderr, argv[0]);
531
532   PCA_Flags = (PCA_CAN_SET_DEFFNM | (MASTER(cr) ? 0 : PCA_QUIET));
533   
534   /* Comment this in to do fexist calls only on master
535    * works not with rerun or tables at the moment
536    * also comment out the version of init_forcerec in md.c 
537    * with NULL instead of opt2fn
538    */
539   /*
540      if (!MASTER(cr))
541      {
542      PCA_Flags |= PCA_NOT_READ_NODE;
543      }
544      */
545
546   parse_common_args(&argc,argv,PCA_Flags, NFILE,fnm,asize(pa),pa,
547                     asize(desc),desc,0,NULL, &oenv);
548
549
550
551   /* we set these early because they might be used in init_multisystem() 
552      Note that there is the potential for npme>nnodes until the number of
553      threads is set later on, if there's thread parallelization. That shouldn't
554      lead to problems. */ 
555   dd_node_order = nenum(ddno_opt);
556   cr->npmenodes = npme;
557
558 #ifndef GMX_THREAD_MPI
559   nthreads=1;
560 #endif
561
562   /* now check the -multi and -multidir option */
563   if (opt2bSet("-multidir", NFILE, fnm))
564   {
565       int i;
566       if (nmultisim > 0)
567       {
568           gmx_fatal(FARGS, "mdrun -multi and -multidir options are mutually exclusive.");
569       }
570       nmultisim = opt2fns(&multidir, "-multidir", NFILE, fnm);
571   }
572
573
574   if (repl_ex_nst != 0 && nmultisim < 2)
575       gmx_fatal(FARGS,"Need at least two replicas for replica exchange (option -multi)");
576
577   if (repl_ex_nex < 0)
578       gmx_fatal(FARGS,"Replica exchange number of exchanges needs to be positive");
579
580   if (nmultisim > 1) {
581 #ifndef GMX_THREAD_MPI
582     gmx_bool bParFn = (multidir == NULL);
583     init_multisystem(cr, nmultisim, multidir, NFILE, fnm, bParFn);
584 #else
585     gmx_fatal(FARGS,"mdrun -multi is not supported with the thread library.Please compile GROMACS with MPI support");
586 #endif
587   }
588
589   bAddPart = !bAppendFiles;
590
591   /* Check if there is ANY checkpoint file available */ 
592   sim_part    = 1;
593   sim_part_fn = sim_part;
594   if (opt2bSet("-cpi",NFILE,fnm))
595   {
596       if (bSepPot && bAppendFiles)
597       {
598           gmx_fatal(FARGS,"Output file appending is not supported with -seppot");
599       }
600
601       bAppendFiles =
602                 read_checkpoint_simulation_part(opt2fn_master("-cpi", NFILE,
603                                                               fnm,cr),
604                                                 &sim_part_fn,NULL,cr,
605                                                 bAppendFiles,NFILE,fnm,
606                                                 part_suffix,&bAddPart);
607       if (sim_part_fn==0 && MASTER(cr))
608       {
609           fprintf(stdout,"No previous checkpoint file present, assuming this is a new run.\n");
610       }
611       else
612       {
613           sim_part = sim_part_fn + 1;
614       }
615
616       if (MULTISIM(cr) && MASTER(cr))
617       {
618           check_multi_int(stdout,cr->ms,sim_part,"simulation part");
619       }
620   } 
621   else
622   {
623       bAppendFiles = FALSE;
624   }
625
626   if (!bAppendFiles)
627   {
628       sim_part_fn = sim_part;
629   }
630
631   if (bAddPart)
632   {
633       /* Rename all output files (except checkpoint files) */
634       /* create new part name first (zero-filled) */
635       sprintf(suffix,"%s%04d",part_suffix,sim_part_fn);
636
637       add_suffix_to_output_names(fnm,NFILE,suffix);
638       if (MASTER(cr))
639       {
640           fprintf(stdout,"Checkpoint file is from part %d, new output files will be suffixed '%s'.\n",sim_part-1,suffix);
641       }
642   }
643
644   Flags = opt2bSet("-rerun",NFILE,fnm) ? MD_RERUN : 0;
645   Flags = Flags | (bSepPot       ? MD_SEPPOT       : 0);
646   Flags = Flags | (bIonize       ? MD_IONIZE       : 0);
647   Flags = Flags | (bPartDec      ? MD_PARTDEC      : 0);
648   Flags = Flags | (bDDBondCheck  ? MD_DDBONDCHECK  : 0);
649   Flags = Flags | (bDDBondComm   ? MD_DDBONDCOMM   : 0);
650   Flags = Flags | (bConfout      ? MD_CONFOUT      : 0);
651   Flags = Flags | (bRerunVSite   ? MD_RERUN_VSITE  : 0);
652   Flags = Flags | (bReproducible ? MD_REPRODUCIBLE : 0);
653   Flags = Flags | (bAppendFiles  ? MD_APPENDFILES  : 0); 
654   Flags = Flags | (opt2parg_bSet("-append", asize(pa),pa) ? MD_APPENDFILESSET : 0); 
655   Flags = Flags | (bKeepAndNumCPT ? MD_KEEPANDNUMCPT : 0); 
656   Flags = Flags | (sim_part>1    ? MD_STARTFROMCPT : 0); 
657   Flags = Flags | (bResetCountersHalfWay ? MD_RESETCOUNTERSHALFWAY : 0);
658
659
660   /* We postpone opening the log file if we are appending, so we can 
661      first truncate the old log file and append to the correct position 
662      there instead.  */
663   if ((MASTER(cr) || bSepPot) && !bAppendFiles) 
664   {
665       gmx_log_open(ftp2fn(efLOG,NFILE,fnm),cr,!bSepPot,Flags,&fplog);
666       CopyRight(fplog,argv[0]);
667       please_cite(fplog,"Hess2008b");
668       please_cite(fplog,"Spoel2005a");
669       please_cite(fplog,"Lindahl2001a");
670       please_cite(fplog,"Berendsen95a");
671   }
672   else if (!MASTER(cr) && bSepPot)
673   {
674       gmx_log_open(ftp2fn(efLOG,NFILE,fnm),cr,!bSepPot,Flags,&fplog);
675   }
676   else
677   {
678       fplog = NULL;
679   }
680
681   ddxyz[XX] = (int)(realddxyz[XX] + 0.5);
682   ddxyz[YY] = (int)(realddxyz[YY] + 0.5);
683   ddxyz[ZZ] = (int)(realddxyz[ZZ] + 0.5);
684
685   rc = mdrunner(nthreads, fplog,cr,NFILE,fnm,oenv,bVerbose,bCompact,
686                 nstglobalcomm, ddxyz,dd_node_order,rdd,rconstr,
687                 dddlb_opt[0],dlb_scale,ddcsx,ddcsy,ddcsz,
688                 nstepout,resetstep,nmultisim,repl_ex_nst,repl_ex_nex,repl_ex_seed,
689                 pforce, cpt_period,max_hours,deviceOptions,Flags);
690
691   gmx_finalize_par();
692
693   if (MULTIMASTER(cr)) {
694       thanx(stderr);
695   }
696
697   /* Log file has to be closed in mdrunner if we are appending to it 
698      (fplog not set here) */
699   if (MASTER(cr) && !bAppendFiles) 
700   {
701       gmx_log_close(fplog);
702   }
703
704   return rc;
705 }
706