Make (most) HTML links to file formats work again
[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,2015, 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 "gmxpre.h"
38
39 #include "config.h"
40
41 #include <stdio.h>
42 #include <string.h>
43
44 #include "gromacs/commandline/pargs.h"
45 #include "gromacs/fileio/filenm.h"
46 #include "gromacs/legacyheaders/checkpoint.h"
47 #include "gromacs/legacyheaders/copyrite.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 #include "gromacs/utility/fatalerror.h"
56
57 #include "mdrun_main.h"
58
59 static bool is_multisim_option_set(int argc, const char *const argv[])
60 {
61     for (int i = 0; i < argc; ++i)
62     {
63         if (strcmp(argv[i], "-multi") == 0 || strcmp(argv[i], "-multidir") == 0)
64         {
65             return true;
66         }
67     }
68     return false;
69 }
70
71 int gmx_mdrun(int argc, char *argv[])
72 {
73     const char   *desc[] = {
74         "[THISMODULE] is the main computational chemistry engine",
75         "within GROMACS. Obviously, it performs Molecular Dynamics simulations,",
76         "but it can also perform Stochastic Dynamics, Energy Minimization,",
77         "test particle insertion or (re)calculation of energies.",
78         "Normal mode analysis is another option. In this case [TT]mdrun[tt]",
79         "builds a Hessian matrix from single conformation.",
80         "For usual Normal Modes-like calculations, make sure that",
81         "the structure provided is properly energy-minimized.",
82         "The generated matrix can be diagonalized by [gmx-nmeig].[PAR]",
83         "The [TT]mdrun[tt] program reads the run input file ([TT]-s[tt])",
84         "and distributes the topology over ranks if needed.",
85         "[TT]mdrun[tt] produces at least four output files.",
86         "A single log file ([TT]-g[tt]) is written.",
87         "The trajectory file ([TT]-o[tt]), contains coordinates, velocities and",
88         "optionally forces.",
89         "The structure file ([TT]-c[tt]) contains the coordinates and",
90         "velocities of the last step.",
91         "The energy file ([TT]-e[tt]) contains energies, the temperature,",
92         "pressure, etc, a lot of these things are also printed in the log file.",
93         "Optionally coordinates can be written to a compressed trajectory file",
94         "([TT]-x[tt]).[PAR]",
95         "The option [TT]-dhdl[tt] is only used when free energy calculation is",
96         "turned on.[PAR]",
97         "Running mdrun efficiently in parallel is a complex topic topic,",
98         "many aspects of which are covered in the online User Guide. You",
99         "should look there for practical advice on using many of the options",
100         "available in mdrun.[PAR]",
101         "ED (essential dynamics) sampling and/or additional flooding potentials",
102         "are switched on by using the [TT]-ei[tt] flag followed by an [REF].edi[ref]",
103         "file. The [REF].edi[ref] file can be produced with the [TT]make_edi[tt] tool",
104         "or by using options in the essdyn menu of the WHAT IF program.",
105         "[TT]mdrun[tt] produces a [REF].xvg[ref] output file that",
106         "contains projections of positions, velocities and forces onto selected",
107         "eigenvectors.[PAR]",
108         "When user-defined potential functions have been selected in the",
109         "[REF].mdp[ref] file the [TT]-table[tt] option is used to pass [TT]mdrun[tt]",
110         "a formatted table with potential functions. The file is read from",
111         "either the current directory or from the [TT]GMXLIB[tt] directory.",
112         "A number of pre-formatted tables are presented in the [TT]GMXLIB[tt] dir,",
113         "for 6-8, 6-9, 6-10, 6-11, 6-12 Lennard-Jones potentials with",
114         "normal Coulomb.",
115         "When pair interactions are present, a separate table for pair interaction",
116         "functions is read using the [TT]-tablep[tt] option.[PAR]",
117         "When tabulated bonded functions are present in the topology,",
118         "interaction functions are read using the [TT]-tableb[tt] option.",
119         "For each different tabulated interaction type the table file name is",
120         "modified in a different way: before the file extension an underscore is",
121         "appended, then a 'b' for bonds, an 'a' for angles or a 'd' for dihedrals",
122         "and finally the table number of the interaction type.[PAR]",
123         "The options [TT]-px[tt] and [TT]-pf[tt] are used for writing pull COM",
124         "coordinates and forces when pulling is selected",
125         "in the [REF].mdp[ref] file.[PAR]",
126         "Finally some experimental algorithms can be tested when the",
127         "appropriate options have been given. Currently under",
128         "investigation are: polarizability.",
129         "[PAR]",
130         "The option [TT]-membed[tt] does what used to be g_membed, i.e. embed",
131         "a protein into a membrane. The data file should contain the options",
132         "that where passed to g_membed before. The [TT]-mn[tt] and [TT]-mp[tt]",
133         "both apply to this as well.",
134         "[PAR]",
135         "The option [TT]-pforce[tt] is useful when you suspect a simulation",
136         "crashes due to too large forces. With this option coordinates and",
137         "forces of atoms with a force larger than a certain value will",
138         "be printed to stderr.",
139         "[PAR]",
140         "Checkpoints containing the complete state of the system are written",
141         "at regular intervals (option [TT]-cpt[tt]) to the file [TT]-cpo[tt],",
142         "unless option [TT]-cpt[tt] is set to -1.",
143         "The previous checkpoint is backed up to [TT]state_prev.cpt[tt] to",
144         "make sure that a recent state of the system is always available,",
145         "even when the simulation is terminated while writing a checkpoint.",
146         "With [TT]-cpnum[tt] all checkpoint files are kept and appended",
147         "with the step number.",
148         "A simulation can be continued by reading the full state from file",
149         "with option [TT]-cpi[tt]. This option is intelligent in the way that",
150         "if no checkpoint file is found, Gromacs just assumes a normal run and",
151         "starts from the first step of the [REF].tpr[ref] file. By default the output",
152         "will be appending to the existing output files. The checkpoint file",
153         "contains checksums of all output files, such that you will never",
154         "loose data when some output files are modified, corrupt or removed.",
155         "There are three scenarios with [TT]-cpi[tt]:[PAR]",
156         "[TT]*[tt] no files with matching names are present: new output files are written[PAR]",
157         "[TT]*[tt] all files are present with names and checksums matching those stored",
158         "in the checkpoint file: files are appended[PAR]",
159         "[TT]*[tt] otherwise no files are modified and a fatal error is generated[PAR]",
160         "With [TT]-noappend[tt] new output files are opened and the simulation",
161         "part number is added to all output file names.",
162         "Note that in all cases the checkpoint file itself is not renamed",
163         "and will be overwritten, unless its name does not match",
164         "the [TT]-cpo[tt] option.",
165         "[PAR]",
166         "With checkpointing the output is appended to previously written",
167         "output files, unless [TT]-noappend[tt] is used or none of the previous",
168         "output files are present (except for the checkpoint file).",
169         "The integrity of the files to be appended is verified using checksums",
170         "which are stored in the checkpoint file. This ensures that output can",
171         "not be mixed up or corrupted due to file appending. When only some",
172         "of the previous output files are present, a fatal error is generated",
173         "and no old output files are modified and no new output files are opened.",
174         "The result with appending will be the same as from a single run.",
175         "The contents will be binary identical, unless you use a different number",
176         "of ranks or dynamic load balancing or the FFT library uses optimizations",
177         "through timing.",
178         "[PAR]",
179         "With option [TT]-maxh[tt] a simulation is terminated and a checkpoint",
180         "file is written at the first neighbor search step where the run time",
181         "exceeds [TT]-maxh[tt]\\*0.99 hours. This option is particularly useful in",
182         "combination with setting [TT]nsteps[tt] to -1 either in the mdp or using the",
183         "similarly named command line option. This results in an infinite run,",
184         "terminated only when the time limit set by [TT]-maxh[tt] is reached (if any)"
185         "or upon receiving a signal."
186         "[PAR]",
187         "When [TT]mdrun[tt] receives a TERM signal, it will set nsteps to the current",
188         "step plus one. When [TT]mdrun[tt] receives an INT signal (e.g. when ctrl+C is",
189         "pressed), it will stop after the next neighbor search step ",
190         "(with nstlist=0 at the next step).",
191         "In both cases all the usual output will be written to file.",
192         "When running with MPI, a signal to one of the [TT]mdrun[tt] ranks",
193         "is sufficient, this signal should not be sent to mpirun or",
194         "the [TT]mdrun[tt] process that is the parent of the others.",
195         "[PAR]",
196         "Interactive molecular dynamics (IMD) can be activated by using at least one",
197         "of the three IMD switches: The [TT]-imdterm[tt] switch allows to terminate the",
198         "simulation from the molecular viewer (e.g. VMD). With [TT]-imdwait[tt],",
199         "[TT]mdrun[tt] pauses whenever no IMD client is connected. Pulling from the",
200         "IMD remote can be turned on by [TT]-imdpull[tt].",
201         "The port [TT]mdrun[tt] listens to can be altered by [TT]-imdport[tt].The",
202         "file pointed to by [TT]-if[tt] contains atom indices and forces if IMD",
203         "pulling is used."
204         "[PAR]",
205         "When [TT]mdrun[tt] is started with MPI, it does not run niced by default."
206     };
207     t_commrec    *cr;
208     t_filenm      fnm[] = {
209         { efTPR, NULL,      NULL,       ffREAD },
210         { efTRN, "-o",      NULL,       ffWRITE },
211         { efCOMPRESSED, "-x", NULL,     ffOPTWR },
212         { efCPT, "-cpi",    NULL,       ffOPTRD | ffALLOW_MISSING },
213         { efCPT, "-cpo",    NULL,       ffOPTWR },
214         { efSTO, "-c",      "confout",  ffWRITE },
215         { efEDR, "-e",      "ener",     ffWRITE },
216         { efLOG, "-g",      "md",       ffWRITE },
217         { efXVG, "-dhdl",   "dhdl",     ffOPTWR },
218         { efXVG, "-field",  "field",    ffOPTWR },
219         { efXVG, "-table",  "table",    ffOPTRD },
220         { efXVG, "-tabletf", "tabletf",    ffOPTRD },
221         { efXVG, "-tablep", "tablep",   ffOPTRD },
222         { efXVG, "-tableb", "table",    ffOPTRD },
223         { efTRX, "-rerun",  "rerun",    ffOPTRD },
224         { efXVG, "-tpi",    "tpi",      ffOPTWR },
225         { efXVG, "-tpid",   "tpidist",  ffOPTWR },
226         { efEDI, "-ei",     "sam",      ffOPTRD },
227         { efXVG, "-eo",     "edsam",    ffOPTWR },
228         { efXVG, "-devout", "deviatie", ffOPTWR },
229         { efXVG, "-runav",  "runaver",  ffOPTWR },
230         { efXVG, "-px",     "pullx",    ffOPTWR },
231         { efXVG, "-pf",     "pullf",    ffOPTWR },
232         { efXVG, "-ro",     "rotation", ffOPTWR },
233         { efLOG, "-ra",     "rotangles", ffOPTWR },
234         { efLOG, "-rs",     "rotslabs", ffOPTWR },
235         { efLOG, "-rt",     "rottorque", ffOPTWR },
236         { efMTX, "-mtx",    "nm",       ffOPTWR },
237         { efNDX, "-dn",     "dipole",   ffOPTWR },
238         { efRND, "-multidir", NULL,      ffOPTRDMULT},
239         { efDAT, "-membed", "membed",   ffOPTRD },
240         { efTOP, "-mp",     "membed",   ffOPTRD },
241         { efNDX, "-mn",     "membed",   ffOPTRD },
242         { efXVG, "-if",     "imdforces", ffOPTWR },
243         { efXVG, "-swap",   "swapions", ffOPTWR }
244     };
245 #define NFILE asize(fnm)
246
247     /* Command line options ! */
248     gmx_bool        bDDBondCheck  = TRUE;
249     gmx_bool        bDDBondComm   = TRUE;
250     gmx_bool        bTunePME      = TRUE;
251     gmx_bool        bVerbose      = FALSE;
252     gmx_bool        bCompact      = TRUE;
253     gmx_bool        bRerunVSite   = FALSE;
254     gmx_bool        bConfout      = TRUE;
255     gmx_bool        bReproducible = FALSE;
256     gmx_bool        bIMDwait      = FALSE;
257     gmx_bool        bIMDterm      = FALSE;
258     gmx_bool        bIMDpull      = FALSE;
259
260     int             npme          = -1;
261     int             nstlist       = 0;
262     int             nmultisim     = 0;
263     int             nstglobalcomm = -1;
264     int             repl_ex_nst   = 0;
265     int             repl_ex_seed  = -1;
266     int             repl_ex_nex   = 0;
267     int             nstepout      = 100;
268     int             resetstep     = -1;
269     gmx_int64_t     nsteps        = -2;   /* the value -2 means that the mdp option will be used */
270     int             imdport       = 8888; /* can be almost anything, 8888 is easy to remember */
271
272     rvec            realddxyz          = {0, 0, 0};
273     const char     *ddno_opt[ddnoNR+1] =
274     { NULL, "interleave", "pp_pme", "cartesian", NULL };
275     const char     *dddlb_opt[] =
276     { NULL, "auto", "no", "yes", NULL };
277     const char     *thread_aff_opt[threadaffNR+1] =
278     { NULL, "auto", "on", "off", NULL };
279     const char     *nbpu_opt[] =
280     { NULL, "auto", "cpu", "gpu", "gpu_cpu", NULL };
281     real            rdd                   = 0.0, rconstr = 0.0, dlb_scale = 0.8, pforce = -1;
282     char           *ddcsx                 = NULL, *ddcsy = NULL, *ddcsz = NULL;
283     real            cpt_period            = 15.0, max_hours = -1;
284     gmx_bool        bAppendFiles          = TRUE;
285     gmx_bool        bKeepAndNumCPT        = FALSE;
286     gmx_bool        bResetCountersHalfWay = FALSE;
287     output_env_t    oenv                  = NULL;
288     const char     *deviceOptions         = "";
289
290     /* Non transparent initialization of a complex gmx_hw_opt_t struct.
291      * But unfortunately we are not allowed to call a function here,
292      * since declarations follow below.
293      */
294     gmx_hw_opt_t    hw_opt = {
295         0, 0, 0, 0, threadaffSEL, 0, 0,
296         { NULL, FALSE, 0, NULL }
297     };
298
299     t_pargs         pa[] = {
300
301         { "-dd",      FALSE, etRVEC, {&realddxyz},
302           "Domain decomposition grid, 0 is optimize" },
303         { "-ddorder", FALSE, etENUM, {ddno_opt},
304           "DD rank order" },
305         { "-npme",    FALSE, etINT, {&npme},
306           "Number of separate ranks to be used for PME, -1 is guess" },
307         { "-nt",      FALSE, etINT, {&hw_opt.nthreads_tot},
308           "Total number of threads to start (0 is guess)" },
309         { "-ntmpi",   FALSE, etINT, {&hw_opt.nthreads_tmpi},
310           "Number of thread-MPI threads to start (0 is guess)" },
311         { "-ntomp",   FALSE, etINT, {&hw_opt.nthreads_omp},
312           "Number of OpenMP threads per MPI rank to start (0 is guess)" },
313         { "-ntomp_pme", FALSE, etINT, {&hw_opt.nthreads_omp_pme},
314           "Number of OpenMP threads per MPI rank to start (0 is -ntomp)" },
315         { "-pin",     FALSE, etENUM, {thread_aff_opt},
316           "Whether mdrun should try to set thread affinities" },
317         { "-pinoffset", FALSE, etINT, {&hw_opt.core_pinning_offset},
318           "The lowest logical core number to which mdrun should pin the first thread" },
319         { "-pinstride", FALSE, etINT, {&hw_opt.core_pinning_stride},
320           "Pinning distance in logical cores for threads, use 0 to minimize the number of threads per physical core" },
321         { "-gpu_id",  FALSE, etSTR, {&hw_opt.gpu_opt.gpu_id},
322           "List of GPU device id-s to use, specifies the per-node PP rank to GPU mapping" },
323         { "-ddcheck", FALSE, etBOOL, {&bDDBondCheck},
324           "Check for all bonded interactions with DD" },
325         { "-ddbondcomm", FALSE, etBOOL, {&bDDBondComm},
326           "HIDDENUse special bonded atom communication when [TT]-rdd[tt] > cut-off" },
327         { "-rdd",     FALSE, etREAL, {&rdd},
328           "The maximum distance for bonded interactions with DD (nm), 0 is determine from initial coordinates" },
329         { "-rcon",    FALSE, etREAL, {&rconstr},
330           "Maximum distance for P-LINCS (nm), 0 is estimate" },
331         { "-dlb",     FALSE, etENUM, {dddlb_opt},
332           "Dynamic load balancing (with DD)" },
333         { "-dds",     FALSE, etREAL, {&dlb_scale},
334           "Fraction in (0,1) by whose reciprocal the initial DD cell size will be increased in order to "
335           "provide a margin in which dynamic load balancing can act while preserving the minimum cell size." },
336         { "-ddcsx",   FALSE, etSTR, {&ddcsx},
337           "HIDDENA string containing a vector of the relative sizes in the x "
338           "direction of the corresponding DD cells. Only effective with static "
339           "load balancing." },
340         { "-ddcsy",   FALSE, etSTR, {&ddcsy},
341           "HIDDENA string containing a vector of the relative sizes in the y "
342           "direction of the corresponding DD cells. Only effective with static "
343           "load balancing." },
344         { "-ddcsz",   FALSE, etSTR, {&ddcsz},
345           "HIDDENA string containing a vector of the relative sizes in the z "
346           "direction of the corresponding DD cells. Only effective with static "
347           "load balancing." },
348         { "-gcom",    FALSE, etINT, {&nstglobalcomm},
349           "Global communication frequency" },
350         { "-nb",      FALSE, etENUM, {&nbpu_opt},
351           "Calculate non-bonded interactions on" },
352         { "-nstlist", FALSE, etINT, {&nstlist},
353           "Set nstlist when using a Verlet buffer tolerance (0 is guess)" },
354         { "-tunepme", FALSE, etBOOL, {&bTunePME},
355           "Optimize PME load between PP/PME ranks or GPU/CPU" },
356         { "-v",       FALSE, etBOOL, {&bVerbose},
357           "Be loud and noisy" },
358         { "-compact", FALSE, etBOOL, {&bCompact},
359           "Write a compact log file" },
360         { "-pforce",  FALSE, etREAL, {&pforce},
361           "Print all forces larger than this (kJ/mol nm)" },
362         { "-reprod",  FALSE, etBOOL, {&bReproducible},
363           "Try to avoid optimizations that affect binary reproducibility" },
364         { "-cpt",     FALSE, etREAL, {&cpt_period},
365           "Checkpoint interval (minutes)" },
366         { "-cpnum",   FALSE, etBOOL, {&bKeepAndNumCPT},
367           "Keep and number checkpoint files" },
368         { "-append",  FALSE, etBOOL, {&bAppendFiles},
369           "Append to previous output files when continuing from checkpoint instead of adding the simulation part number to all file names" },
370         { "-nsteps",  FALSE, etINT64, {&nsteps},
371           "Run this number of steps, overrides .mdp file option (-1 means infinite, -2 means use mdp option, smaller is invalid)" },
372         { "-maxh",   FALSE, etREAL, {&max_hours},
373           "Terminate after 0.99 times this time (hours)" },
374         { "-multi",   FALSE, etINT, {&nmultisim},
375           "Do multiple simulations in parallel" },
376         { "-replex",  FALSE, etINT, {&repl_ex_nst},
377           "Attempt replica exchange periodically with this period (steps)" },
378         { "-nex",  FALSE, etINT, {&repl_ex_nex},
379           "Number of random exchanges to carry out each exchange interval (N^3 is one suggestion).  -nex zero or not specified gives neighbor replica exchange." },
380         { "-reseed",  FALSE, etINT, {&repl_ex_seed},
381           "Seed for replica exchange, -1 is generate a seed" },
382         { "-imdport",    FALSE, etINT, {&imdport},
383           "HIDDENIMD listening port" },
384         { "-imdwait",  FALSE, etBOOL, {&bIMDwait},
385           "HIDDENPause the simulation while no IMD client is connected" },
386         { "-imdterm",  FALSE, etBOOL, {&bIMDterm},
387           "HIDDENAllow termination of the simulation from IMD client" },
388         { "-imdpull",  FALSE, etBOOL, {&bIMDpull},
389           "HIDDENAllow pulling in the simulation from IMD client" },
390         { "-rerunvsite", FALSE, etBOOL, {&bRerunVSite},
391           "HIDDENRecalculate virtual site coordinates with [TT]-rerun[tt]" },
392         { "-confout", FALSE, etBOOL, {&bConfout},
393           "HIDDENWrite the last configuration with [TT]-c[tt] and force checkpointing at the last step" },
394         { "-stepout", FALSE, etINT, {&nstepout},
395           "HIDDENFrequency of writing the remaining wall clock time for the run" },
396         { "-resetstep", FALSE, etINT, {&resetstep},
397           "HIDDENReset cycle counters after these many time steps" },
398         { "-resethway", FALSE, etBOOL, {&bResetCountersHalfWay},
399           "HIDDENReset the cycle counters after half the number of steps or halfway [TT]-maxh[tt]" }
400     };
401     unsigned long   Flags;
402     ivec            ddxyz;
403     int             dd_node_order;
404     gmx_bool        bAddPart;
405     FILE           *fplog, *fpmulti;
406     int             sim_part, sim_part_fn;
407     const char     *part_suffix = ".part";
408     char            suffix[STRLEN];
409     int             rc;
410     char          **multidir = NULL;
411
412     cr = init_commrec();
413
414     unsigned long PCA_Flags = PCA_CAN_SET_DEFFNM;
415     // With -multi or -multidir, the file names are going to get processed
416     // further (or the working directory changed), so we can't check for their
417     // existence during parsing.  It isn't useful to do any completion based on
418     // file system contents, either.
419     if (is_multisim_option_set(argc, argv))
420     {
421         PCA_Flags |= PCA_DISABLE_INPUT_FILE_CHECKING;
422     }
423
424     /* Comment this in to do fexist calls only on master
425      * works not with rerun or tables at the moment
426      * also comment out the version of init_forcerec in md.c
427      * with NULL instead of opt2fn
428      */
429     /*
430        if (!MASTER(cr))
431        {
432        PCA_Flags |= PCA_NOT_READ_NODE;
433        }
434      */
435
436     if (!parse_common_args(&argc, argv, PCA_Flags, NFILE, fnm, asize(pa), pa,
437                            asize(desc), desc, 0, NULL, &oenv))
438     {
439         return 0;
440     }
441
442
443     /* we set these early because they might be used in init_multisystem()
444        Note that there is the potential for npme>nnodes until the number of
445        threads is set later on, if there's thread parallelization. That shouldn't
446        lead to problems. */
447     dd_node_order = nenum(ddno_opt);
448     cr->npmenodes = npme;
449
450     hw_opt.thread_affinity = nenum(thread_aff_opt);
451
452     /* now check the -multi and -multidir option */
453     if (opt2bSet("-multidir", NFILE, fnm))
454     {
455         if (nmultisim > 0)
456         {
457             gmx_fatal(FARGS, "mdrun -multi and -multidir options are mutually exclusive.");
458         }
459         nmultisim = opt2fns(&multidir, "-multidir", NFILE, fnm);
460     }
461
462
463     if (repl_ex_nst != 0 && nmultisim < 2)
464     {
465         gmx_fatal(FARGS, "Need at least two replicas for replica exchange (option -multi)");
466     }
467
468     if (repl_ex_nex < 0)
469     {
470         gmx_fatal(FARGS, "Replica exchange number of exchanges needs to be positive");
471     }
472
473     if (nmultisim > 1)
474     {
475 #ifndef GMX_THREAD_MPI
476         gmx_bool bParFn = (multidir == NULL);
477         init_multisystem(cr, nmultisim, multidir, NFILE, fnm, bParFn);
478 #else
479         gmx_fatal(FARGS, "mdrun -multi is not supported with the thread library. "
480                   "Please compile GROMACS with MPI support");
481 #endif
482     }
483
484     bAddPart = !bAppendFiles;
485
486     /* Check if there is ANY checkpoint file available */
487     sim_part    = 1;
488     sim_part_fn = sim_part;
489     if (opt2bSet("-cpi", NFILE, fnm))
490     {
491         bAppendFiles =
492             read_checkpoint_simulation_part(opt2fn_master("-cpi", NFILE,
493                                                           fnm, cr),
494                                             &sim_part_fn, NULL, cr,
495                                             bAppendFiles, NFILE, fnm,
496                                             part_suffix, &bAddPart);
497         if (sim_part_fn == 0 && MULTIMASTER(cr))
498         {
499             fprintf(stdout, "No previous checkpoint file present, assuming this is a new run.\n");
500         }
501         else
502         {
503             sim_part = sim_part_fn + 1;
504         }
505
506         if (MULTISIM(cr) && MASTER(cr))
507         {
508             if (MULTIMASTER(cr))
509             {
510                 /* Log file is not yet available, so if there's a
511                  * problem we can only write to stderr. */
512                 fpmulti = stderr;
513             }
514             else
515             {
516                 fpmulti = NULL;
517             }
518             check_multi_int(fpmulti, cr->ms, sim_part, "simulation part", TRUE);
519         }
520     }
521     else
522     {
523         bAppendFiles = FALSE;
524     }
525
526     if (!bAppendFiles)
527     {
528         sim_part_fn = sim_part;
529     }
530
531     if (bAddPart)
532     {
533         /* Rename all output files (except checkpoint files) */
534         /* create new part name first (zero-filled) */
535         sprintf(suffix, "%s%04d", part_suffix, sim_part_fn);
536
537         add_suffix_to_output_names(fnm, NFILE, suffix);
538         if (MULTIMASTER(cr))
539         {
540             fprintf(stdout, "Checkpoint file is from part %d, new output files will be suffixed '%s'.\n", sim_part-1, suffix);
541         }
542     }
543
544     Flags = opt2bSet("-rerun", NFILE, fnm) ? MD_RERUN : 0;
545     Flags = Flags | (bDDBondCheck  ? MD_DDBONDCHECK  : 0);
546     Flags = Flags | (bDDBondComm   ? MD_DDBONDCOMM   : 0);
547     Flags = Flags | (bTunePME      ? MD_TUNEPME      : 0);
548     Flags = Flags | (bConfout      ? MD_CONFOUT      : 0);
549     Flags = Flags | (bRerunVSite   ? MD_RERUN_VSITE  : 0);
550     Flags = Flags | (bReproducible ? MD_REPRODUCIBLE : 0);
551     Flags = Flags | (bAppendFiles  ? MD_APPENDFILES  : 0);
552     Flags = Flags | (opt2parg_bSet("-append", asize(pa), pa) ? MD_APPENDFILESSET : 0);
553     Flags = Flags | (bKeepAndNumCPT ? MD_KEEPANDNUMCPT : 0);
554     Flags = Flags | (sim_part > 1    ? MD_STARTFROMCPT : 0);
555     Flags = Flags | (bResetCountersHalfWay ? MD_RESETCOUNTERSHALFWAY : 0);
556     Flags = Flags | (bIMDwait      ? MD_IMDWAIT      : 0);
557     Flags = Flags | (bIMDterm      ? MD_IMDTERM      : 0);
558     Flags = Flags | (bIMDpull      ? MD_IMDPULL      : 0);
559
560     /* We postpone opening the log file if we are appending, so we can
561        first truncate the old log file and append to the correct position
562        there instead.  */
563     if (MASTER(cr) && !bAppendFiles)
564     {
565         gmx_log_open(ftp2fn(efLOG, NFILE, fnm), cr,
566                      Flags & MD_APPENDFILES, &fplog);
567         please_cite(fplog, "Hess2008b");
568         please_cite(fplog, "Spoel2005a");
569         please_cite(fplog, "Lindahl2001a");
570         please_cite(fplog, "Berendsen95a");
571     }
572     else
573     {
574         fplog = NULL;
575     }
576
577     ddxyz[XX] = (int)(realddxyz[XX] + 0.5);
578     ddxyz[YY] = (int)(realddxyz[YY] + 0.5);
579     ddxyz[ZZ] = (int)(realddxyz[ZZ] + 0.5);
580
581     rc = mdrunner(&hw_opt, fplog, cr, NFILE, fnm, oenv, bVerbose, bCompact,
582                   nstglobalcomm, ddxyz, dd_node_order, rdd, rconstr,
583                   dddlb_opt[0], dlb_scale, ddcsx, ddcsy, ddcsz,
584                   nbpu_opt[0], nstlist,
585                   nsteps, nstepout, resetstep,
586                   nmultisim, repl_ex_nst, repl_ex_nex, repl_ex_seed,
587                   pforce, cpt_period, max_hours, deviceOptions, imdport, Flags);
588
589     /* Log file has to be closed in mdrunner if we are appending to it
590        (fplog not set here) */
591     if (MASTER(cr) && !bAppendFiles)
592     {
593         gmx_log_close(fplog);
594     }
595
596     return rc;
597 }