made mdrun file appending default and added checks for presence and suffix of output...
[gromacs/rigid-bodies.git] / src / kernel / mdrun.c
blob03f46db80c0ca63d12cabe9863f4a061ea75b1eb
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
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.
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.
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.
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.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
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"
53 /* afm stuf */
54 #include "pull.h"
56 int main(int argc,char *argv[])
58 const char *desc[] = {
59 "The mdrun program is the main computational chemistry engine",
60 "within GROMACS. Obviously, it performs Molecular Dynamics simulations,",
61 "but it can also perform Stochastic Dynamics, Energy Minimization,",
62 "test particle insertion or (re)calculation of energies.",
63 "Normal mode analysis is another option. In this case mdrun",
64 "builds a Hessian matrix from single conformation.",
65 "For usual Normal Modes-like calculations, make sure that",
66 "the structure provided is properly energy-minimized.",
67 "The generated matrix can be diagonalized by g_nmeig.[PAR]",
68 "The mdrun program reads the run input file ([TT]-s[tt])",
69 "and distributes the topology over nodes if needed.",
70 "mdrun produces at least four output files.",
71 "A single log file ([TT]-g[tt]) is written, unless the option",
72 "[TT]-seppot[tt] is used, in which case each node writes a log file.",
73 "The trajectory file ([TT]-o[tt]), contains coordinates, velocities and",
74 "optionally forces.",
75 "The structure file ([TT]-c[tt]) contains the coordinates and",
76 "velocities of the last step.",
77 "The energy file ([TT]-e[tt]) contains energies, the temperature,",
78 "pressure, etc, a lot of these things are also printed in the log file.",
79 "Optionally coordinates can be written to a compressed trajectory file",
80 "([TT]-x[tt]).[PAR]",
81 "The option [TT]-dhdl[tt] is only used when free energy calculation is",
82 "turned on.[PAR]",
83 "When mdrun is started using MPI with more than 1 node, parallelization",
84 "is used. By default domain decomposition is used, unless the [TT]-pd[tt]",
85 "option is set, which selects particle decomposition.[PAR]",
86 "With domain decomposition, the spatial decomposition can be set",
87 "with option [TT]-dd[tt]. By default mdrun selects a good decomposition.",
88 "The user only needs to change this when the system is very inhomogeneous.",
89 "Dynamic load balancing is set with the option [TT]-dlb[tt],",
90 "which can give a significant performance improvement,",
91 "especially for inhomogeneous systems. The only disadvantage of",
92 "dynamic load balancing is that runs are no longer binary reproducible,",
93 "but in most cases this is not important.",
94 "By default the dynamic load balancing is automatically turned on",
95 "when the measured performance loss due to load imbalance is 5% or more.",
96 "At low parallelization these are the only important options",
97 "for domain decomposition.",
98 "At high parallelization the options in the next two sections",
99 "could be important for increasing the performace.",
100 "[PAR]",
101 "When PME is used with domain decomposition, separate nodes can",
102 "be assigned to do only the PME mesh calculation;",
103 "this is computationally more efficient starting at about 12 nodes.",
104 "The number of PME nodes is set with option [TT]-npme[tt],",
105 "this can not be more than half of the nodes.",
106 "By default mdrun makes a guess for the number of PME",
107 "nodes when the number of nodes is larger than 11 or performance wise",
108 "not compatible with the PME grid x dimension.",
109 "But the user should optimize npme. Performance statistics on this issue",
110 "are written at the end of the log file.",
111 "For good load balancing at high parallelization, the PME grid x and y",
112 "dimensions should be divisible by the number of PME nodes",
113 "(the simulation will run correctly also when this is not the case).",
114 "[PAR]",
115 "This section lists all options that affect the domain decomposition.",
116 "[BR]",
117 "Option [TT]-rdd[tt] can be used to set the required maximum distance",
118 "for inter charge-group bonded interactions.",
119 "Communication for two-body bonded interactions below the non-bonded",
120 "cut-off distance always comes for free with the non-bonded communication.",
121 "Atoms beyond the non-bonded cut-off are only communicated when they have",
122 "missing bonded interactions; this means that the extra cost is minor",
123 "and nearly indepedent of the value of [TT]-rdd[tt].",
124 "With dynamic load balancing option [TT]-rdd[tt] also sets",
125 "the lower limit for the domain decomposition cell sizes.",
126 "By default [TT]-rdd[tt] is determined by mdrun based on",
127 "the initial coordinates. The chosen value will be a balance",
128 "between interaction range and communication cost.",
129 "[BR]",
130 "When inter charge-group bonded interactions are beyond",
131 "the bonded cut-off distance, mdrun terminates with an error message.",
132 "For pair interactions and tabulated bonds",
133 "that do not generate exclusions, this check can be turned off",
134 "with the option [TT]-noddcheck[tt].",
135 "[BR]",
136 "When constraints are present, option [TT]-rcon[tt] influences",
137 "the cell size limit as well.",
138 "Atoms connected by NC constraints, where NC is the LINCS order plus 1,",
139 "should not be beyond the smallest cell size. A error message is",
140 "generated when this happens and the user should change the decomposition",
141 "or decrease the LINCS order and increase the number of LINCS iterations.",
142 "By default mdrun estimates the minimum cell size required for P-LINCS",
143 "in a conservative fashion. For high parallelization it can be useful",
144 "to set the distance required for P-LINCS with the option [TT]-rcon[tt].",
145 "[BR]",
146 "The [TT]-dds[tt] option sets the minimum allowed x, y and/or z scaling",
147 "of the cells with dynamic load balancing. mdrun will ensure that",
148 "the cells can scale down by at least this factor. This option is used",
149 "for the automated spatial decomposition (when not using [TT]-dd[tt])",
150 "as well as for determining the number of grid pulses, which in turn",
151 "sets the minimum allowed cell size. Under certain circumstances",
152 "the value of [TT]-dds[tt] might need to be adjusted to account for",
153 "high or low spatial inhomogeneity of the system.",
154 "[PAR]",
155 "The option [TT]-gcom[tt] can be used to only do global communication",
156 "every n steps."
157 "This can improve performance for highly parallel simulations",
158 "where this global communication step becomes the bottleneck.",
159 "For a global thermostat and/or barostat the temperature",
160 "and/or pressure will also only be updated every -gcom steps.",
161 "By default it is set to the minimum of nstcalcenergy and nstlist.[PAR]",
162 "With [TT]-rerun[tt] an input trajectory can be given for which ",
163 "forces and energies will be (re)calculated. Neighbor searching will be",
164 "performed for every frame, unless [TT]nstlist[tt] is zero",
165 "(see the [TT].mdp[tt] file).[PAR]",
166 "ED (essential dynamics) sampling is switched on by using the [TT]-ei[tt]",
167 "flag followed by an [TT].edi[tt] file.",
168 "The [TT].edi[tt] file can be produced using options in the essdyn",
169 "menu of the WHAT IF program. mdrun produces a [TT].edo[tt] file that",
170 "contains projections of positions, velocities and forces onto selected",
171 "eigenvectors.[PAR]",
172 "When user-defined potential functions have been selected in the",
173 "[TT].mdp[tt] file the [TT]-table[tt] option is used to pass mdrun",
174 "a formatted table with potential functions. The file is read from",
175 "either the current directory or from the GMXLIB directory.",
176 "A number of pre-formatted tables are presented in the GMXLIB dir,",
177 "for 6-8, 6-9, 6-10, 6-11, 6-12 Lennard Jones potentials with",
178 "normal Coulomb.",
179 "When pair interactions are present a separate table for pair interaction",
180 "functions is read using the [TT]-tablep[tt] option.[PAR]",
181 "When tabulated bonded functions are present in the topology,",
182 "interaction functions are read using the [TT]-tableb[tt] option.",
183 "For each different tabulated interaction type the table file name is",
184 "modified in a different way: before the file extension an underscore is",
185 "appended, then a b for bonds, an a for angles or a d for dihedrals",
186 "and finally the table number of the interaction type.[PAR]",
187 "The options [TT]-px[tt] and [TT]-pf[tt] are used for writing pull COM",
188 "coordinates and forces when pulling is selected",
189 "in the [TT].mdp[tt] file.[PAR]",
190 "With [TT]-multi[tt] multiple systems are simulated in parallel.",
191 "As many input files are required as the number of systems.",
192 "The system number is appended to the run input and each output filename,",
193 "for instance topol.tpr becomes topol0.tpr, topol1.tpr etc.",
194 "The number of nodes per system is the total number of nodes",
195 "divided by the number of systems.",
196 "One use of this option is for NMR refinement: when distance",
197 "or orientation restraints are present these can be ensemble averaged",
198 "over all the systems.[PAR]",
199 "With [TT]-replex[tt] replica exchange is attempted every given number",
200 "of steps. The number of replicas is set with the [TT]-multi[tt] option,",
201 "see above.",
202 "All run input files should use a different coupling temperature,",
203 "the order of the files is not important. The random seed is set with",
204 "[TT]-reseed[tt]. The velocities are scaled and neighbor searching",
205 "is performed after every exchange.[PAR]",
206 "Finally some experimental algorithms can be tested when the",
207 "appropriate options have been given. Currently under",
208 "investigation are: polarizability, and X-Ray bombardments.",
209 "[PAR]",
210 "The option [TT]-pforce[tt] is useful when you suspect a simulation",
211 "crashes due to too large forces. With this option coordinates and",
212 "forces of atoms with a force larger than a certain value will",
213 "be printed to stderr.",
214 "[PAR]",
215 "Checkpoints containing the complete state of the system are written",
216 "at regular intervals (option [TT]-cpt[tt]) to the file [TT]-cpo[tt],",
217 "unless option [TT]-cpt[tt] is set to -1.",
218 "The previous checkpoint is backed up to [TT]state_prev.cpt[tt] to",
219 "make sure that a recent state of the system is always available,",
220 "even when the simulation is terminated while writing a checkpoint.",
221 "A simulation can be continued by reading the full state from file",
222 "with option [TT]-cpi[tt]. This option is intelligent in the way that",
223 "if no checkpoint file is found, Gromacs just assumes a normal run and",
224 "starts from the first step of the tpr file.",
225 "The simulation part number is added to all output files,",
226 "unless [TT]-append[tt] or [TT]-noaddpart[tt] are set.",
227 "[PAR]",
228 "With checkpointing the output is appended to previously written",
229 "output files, unless [TT]-noappend[tt] is used or none of the previous",
230 "output files are present (except for the checkpoint file).",
231 "The integrity of the files to be appended is verified using checksums",
232 "which are stored in the checkpoint file. This ensures that output can",
233 "not be mixed up or corrupted due to file appending. When only some",
234 "of the previous output files are present, a fatal error is generated",
235 "and no old output files are modified and no new output files are opened.",
236 "The result with appending will be the same as from a single run.",
237 "The contents will be binary identical, unless you use a different number",
238 "of nodes or dynamic load balancing or the FFT library uses optimizations",
239 "through timing."
240 "[PAR]",
241 "With option [TT]-maxh[tt] a simulation is terminated and a checkpoint",
242 "file is written at the first neighbor search step where the run time",
243 "exceeds [TT]-maxh[tt]*0.99 hours.",
244 "[PAR]",
245 "When mdrun receives a TERM signal, it will set nsteps to the current",
246 "step plus one. When mdrun receives an INT signal (e.g. when ctrl+C is",
247 "pressed), it will stop after the next neighbor search step ",
248 "(with nstlist=0 at the next step).",
249 "In both cases all the usual output will be written to file.",
250 "When running with MPI, a signal to one of the mdrun processes",
251 "is sufficient, this signal should not be sent to mpirun or",
252 "the mdrun process that is the parent of the others.",
253 "[PAR]",
254 "When mdrun is started with MPI, it does not run niced by default."
256 t_commrec *cr;
257 t_filenm fnm[] = {
258 { efTPX, NULL, NULL, ffREAD },
259 { efTRN, "-o", NULL, ffWRITE },
260 { efXTC, "-x", NULL, ffOPTWR },
261 { efCPT, "-cpi", NULL, ffOPTRD },
262 { efCPT, "-cpo", NULL, ffOPTWR },
263 { efSTO, "-c", "confout", ffWRITE },
264 { efEDR, "-e", "ener", ffWRITE },
265 { efLOG, "-g", "md", ffWRITE },
266 { efXVG, "-dhdl", "dhdl", ffOPTWR },
267 { efXVG, "-field", "field", ffOPTWR },
268 { efXVG, "-table", "table", ffOPTRD },
269 { efXVG, "-tablep", "tablep", ffOPTRD },
270 { efXVG, "-tableb", "table", ffOPTRD },
271 { efTRX, "-rerun", "rerun", ffOPTRD },
272 { efXVG, "-tpi", "tpi", ffOPTWR },
273 { efXVG, "-tpid", "tpidist", ffOPTWR },
274 { efEDI, "-ei", "sam", ffOPTRD },
275 { efEDO, "-eo", "sam", ffOPTWR },
276 { efGCT, "-j", "wham", ffOPTRD },
277 { efGCT, "-jo", "bam", ffOPTWR },
278 { efXVG, "-ffout", "gct", ffOPTWR },
279 { efXVG, "-devout", "deviatie", ffOPTWR },
280 { efXVG, "-runav", "runaver", ffOPTWR },
281 { efXVG, "-px", "pullx", ffOPTWR },
282 { efXVG, "-pf", "pullf", ffOPTWR },
283 { efMTX, "-mtx", "nm", ffOPTWR },
284 { efNDX, "-dn", "dipole", ffOPTWR }
286 #define NFILE asize(fnm)
288 /* Command line options ! */
289 bool bCart = FALSE;
290 bool bPPPME = FALSE;
291 bool bPartDec = FALSE;
292 bool bDDBondCheck = TRUE;
293 bool bDDBondComm = TRUE;
294 bool bVerbose = FALSE;
295 bool bCompact = TRUE;
296 bool bSepPot = FALSE;
297 bool bRerunVSite = FALSE;
298 bool bIonize = FALSE;
299 bool bConfout = TRUE;
300 bool bReproducible = FALSE;
302 int npme=-1;
303 int nmultisim=0;
304 int nstglobalcomm=-1;
305 int repl_ex_nst=0;
306 int repl_ex_seed=-1;
307 int nstepout=100;
308 int nthreads=1;
309 int resetstep=-1;
311 rvec realddxyz={0,0,0};
312 const char *ddno_opt[ddnoNR+1] =
313 { NULL, "interleave", "pp_pme", "cartesian", NULL };
314 const char *dddlb_opt[] =
315 { NULL, "auto", "no", "yes", NULL };
316 real rdd=0.0,rconstr=0.0,dlb_scale=0.8,pforce=-1;
317 char *ddcsx=NULL,*ddcsy=NULL,*ddcsz=NULL;
318 real cpt_period=15.0,max_hours=-1;
319 bool bAppendFiles=TRUE,bAddPart=TRUE;
320 bool bResetCountersHalfWay=FALSE;
321 output_env_t oenv=NULL;
323 t_pargs pa[] = {
324 { "-pd", FALSE, etBOOL,{&bPartDec},
325 "Use particle decompostion" },
326 { "-dd", FALSE, etRVEC,{&realddxyz},
327 "Domain decomposition grid, 0 is optimize" },
328 { "-nt", FALSE, etINT, {&nthreads},
329 "Number of threads to start on each node" },
330 { "-npme", FALSE, etINT, {&npme},
331 "Number of separate nodes to be used for PME, -1 is guess" },
332 { "-ddorder", FALSE, etENUM, {ddno_opt},
333 "DD node order" },
334 { "-ddcheck", FALSE, etBOOL, {&bDDBondCheck},
335 "Check for all bonded interactions with DD" },
336 { "-ddbondcomm", FALSE, etBOOL, {&bDDBondComm},
337 "HIDDENUse special bonded atom communication when -rdd > cut-off" },
338 { "-rdd", FALSE, etREAL, {&rdd},
339 "The maximum distance for bonded interactions with DD (nm), 0 is determine from initial coordinates" },
340 { "-rcon", FALSE, etREAL, {&rconstr},
341 "Maximum distance for P-LINCS (nm), 0 is estimate" },
342 { "-dlb", FALSE, etENUM, {dddlb_opt},
343 "Dynamic load balancing (with DD)" },
344 { "-dds", FALSE, etREAL, {&dlb_scale},
345 "Minimum allowed dlb scaling of the DD cell size" },
346 { "-ddcsx", FALSE, etSTR, {&ddcsx},
347 "HIDDENThe DD cell sizes in x" },
348 { "-ddcsy", FALSE, etSTR, {&ddcsy},
349 "HIDDENThe DD cell sizes in y" },
350 { "-ddcsz", FALSE, etSTR, {&ddcsz},
351 "HIDDENThe DD cell sizes in z" },
352 { "-gcom", FALSE, etINT,{&nstglobalcomm},
353 "Global communication frequency" },
354 { "-v", FALSE, etBOOL,{&bVerbose},
355 "Be loud and noisy" },
356 { "-compact", FALSE, etBOOL,{&bCompact},
357 "Write a compact log file" },
358 { "-seppot", FALSE, etBOOL, {&bSepPot},
359 "Write separate V and dVdl terms for each interaction type and node to the log file(s)" },
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 { "-append", FALSE, etBOOL, {&bAppendFiles},
367 "Append to previous output files when continuing from checkpoint" },
368 { "-addpart", FALSE, etBOOL, {&bAddPart},
369 "Add the simulation part number to all output files when continuing from checkpoint" },
370 { "-maxh", FALSE, etREAL, {&max_hours},
371 "Terminate after 0.99 times this time (hours)" },
372 { "-multi", FALSE, etINT,{&nmultisim},
373 "Do multiple simulations in parallel" },
374 { "-replex", FALSE, etINT, {&repl_ex_nst},
375 "Attempt replica exchange every # steps" },
376 { "-reseed", FALSE, etINT, {&repl_ex_seed},
377 "Seed for replica exchange, -1 is generate a seed" },
378 { "-rerunvsite", FALSE, etBOOL, {&bRerunVSite},
379 "HIDDENRecalculate virtual site coordinates with -rerun" },
380 { "-ionize", FALSE, etBOOL,{&bIonize},
381 "Do a simulation including the effect of an X-Ray bombardment on your system" },
382 { "-confout", FALSE, etBOOL, {&bConfout},
383 "HIDDENWrite the last configuration with -c and force checkpointing at the last step" },
384 { "-stepout", FALSE, etINT, {&nstepout},
385 "HIDDENFrequency of writing the remaining runtime" },
386 { "-resetstep", FALSE, etINT, {&resetstep},
387 "HIDDENReset cycle counters after these many time steps" },
388 { "-resethway", FALSE, etBOOL, {&bResetCountersHalfWay},
389 "HIDDENReset the cycle counters after half the number of steps or halfway -maxh" }
391 gmx_edsam_t ed;
392 unsigned long Flags, PCA_Flags;
393 ivec ddxyz;
394 int dd_node_order;
395 bool HaveCheckpoint;
396 FILE *fplog,*fptest;
397 int sim_part,sim_part_fn;
398 const char *part_suffix=".part";
399 char suffix[STRLEN];
400 int rc;
402 cr = init_par(&argc,&argv);
403 cr->nthreads = nthreads;
405 PCA_Flags = (PCA_KEEP_ARGS | PCA_NOEXIT_ON_ARGS | PCA_CAN_SET_DEFFNM
406 | (MASTER(cr) ? 0 : PCA_QUIET));
409 /* Comment this in to do fexist calls only on master
410 * works not with rerun or tables at the moment
411 * also comment out the version of init_forcerec in md.c
412 * with NULL instead of opt2fn
415 if (!MASTER(cr))
417 PCA_Flags |= PCA_NOT_READ_NODE;
421 parse_common_args(&argc,argv,PCA_Flags, NFILE,fnm,asize(pa),pa,
422 asize(desc),desc,0,NULL, &oenv);
425 dd_node_order = nenum(ddno_opt);
426 cr->npmenodes = npme;
428 if (repl_ex_nst != 0 && nmultisim < 2)
429 gmx_fatal(FARGS,"Need at least two replicas for replica exchange (option -multi)");
431 if (nmultisim > 1) {
432 #ifndef GMX_THREADS
433 init_multisystem(cr,nmultisim,NFILE,fnm,TRUE);
434 #else
435 gmx_fatal(FARGS,"mdrun -multi is not supported with the thread library.Please compile GROMACS with MPI support");
436 #endif
439 /* Check if there is ANY checkpoint file available */
440 sim_part = 1;
441 sim_part_fn = sim_part;
442 if (opt2bSet("-cpi",NFILE,fnm))
444 bAppendFiles =
445 read_checkpoint_simulation_part(opt2fn_master("-cpi",NFILE,fnm,cr),
446 &sim_part_fn,NULL,cr,
447 bAppendFiles,
448 part_suffix,&bAddPart);
449 if (sim_part_fn==0 && MASTER(cr))
451 fprintf(stdout,"No previous checkpoint file present, assuming this is a new run.\n");
453 else
455 sim_part = sim_part_fn + 1;
458 else
460 bAppendFiles = FALSE;
463 if (!bAppendFiles)
465 sim_part_fn = sim_part;
468 if (bAddPart && sim_part_fn > 1)
470 /* This is a continuation run, rename trajectory output files
471 (except checkpoint files) */
472 /* create new part name first (zero-filled) */
473 sprintf(suffix,"%s%04d",part_suffix,sim_part_fn);
475 add_suffix_to_output_names(fnm,NFILE,suffix);
476 fprintf(stdout,"Checkpoint file is from part %d, new output files will be suffixed '%s'.\n",sim_part-1,suffix);
479 Flags = opt2bSet("-rerun",NFILE,fnm) ? MD_RERUN : 0;
480 Flags = Flags | (bSepPot ? MD_SEPPOT : 0);
481 Flags = Flags | (bIonize ? MD_IONIZE : 0);
482 Flags = Flags | (bPartDec ? MD_PARTDEC : 0);
483 Flags = Flags | (bDDBondCheck ? MD_DDBONDCHECK : 0);
484 Flags = Flags | (bDDBondComm ? MD_DDBONDCOMM : 0);
485 Flags = Flags | (bConfout ? MD_CONFOUT : 0);
486 Flags = Flags | (bRerunVSite ? MD_RERUN_VSITE : 0);
487 Flags = Flags | (bReproducible ? MD_REPRODUCIBLE : 0);
488 Flags = Flags | (bAppendFiles ? MD_APPENDFILES : 0);
489 Flags = Flags | (sim_part>1 ? MD_STARTFROMCPT : 0);
490 Flags = Flags | (bResetCountersHalfWay ? MD_RESETCOUNTERSHALFWAY : 0);
493 /* We postpone opening the log file if we are appending, so we can
494 first truncate the old log file and append to the correct position
495 there instead. */
496 if ((MASTER(cr) || bSepPot) && !bAppendFiles)
498 gmx_log_open(ftp2fn(efLOG,NFILE,fnm),cr,!bSepPot,Flags,&fplog);
499 CopyRight(fplog,argv[0]);
500 please_cite(fplog,"Hess2008b");
501 please_cite(fplog,"Spoel2005a");
502 please_cite(fplog,"Lindahl2001a");
503 please_cite(fplog,"Berendsen95a");
505 else
507 fplog = NULL;
510 #if 0
511 /* this is now done in mdrunner: */
512 /* Essential dynamics */
513 if (opt2bSet("-ei",NFILE,fnm)) {
514 /* Open input and output files, allocate space for ED data structure */
515 ed = ed_open(NFILE,fnm,cr);
516 } else
517 ed=NULL;
518 #endif
520 ddxyz[XX] = (int)(realddxyz[XX] + 0.5);
521 ddxyz[YY] = (int)(realddxyz[YY] + 0.5);
522 ddxyz[ZZ] = (int)(realddxyz[ZZ] + 0.5);
524 /* even if nthreads = 1, we still call this one */
525 rc = mdrunner_threads(nthreads,
526 fplog,cr,NFILE,fnm,oenv,bVerbose,bCompact,nstglobalcomm,
527 ddxyz,dd_node_order,rdd,rconstr,
528 dddlb_opt[0],dlb_scale,ddcsx,ddcsy,ddcsz,
529 nstepout,resetstep,nmultisim,repl_ex_nst,repl_ex_seed,pforce,
530 cpt_period,max_hours,Flags);
532 if (gmx_parallel_env_initialized())
533 gmx_finalize();
535 if (MULTIMASTER(cr)) {
536 thanx(stderr);
539 /* Log file has to be closed in mdrunner if we are appending to it
540 (fplog not set here) */
541 if (MASTER(cr) && !bAppendFiles)
543 gmx_log_close(fplog);
546 return rc;