Merge branch release-4-6 into master
[alexxy/gromacs.git] / src / programs / mdrun / md.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include "typedefs.h"
40 #include "smalloc.h"
41 #include "sysstuff.h"
42 #include "vec.h"
43 #include "statutil.h"
44 #include "vcm.h"
45 #include "mdebin.h"
46 #include "nrnb.h"
47 #include "calcmu.h"
48 #include "index.h"
49 #include "vsite.h"
50 #include "update.h"
51 #include "ns.h"
52 #include "trnio.h"
53 #include "xtcio.h"
54 #include "mdrun.h"
55 #include "md_support.h"
56 #include "md_logging.h"
57 #include "confio.h"
58 #include "network.h"
59 #include "pull.h"
60 #include "xvgr.h"
61 #include "physics.h"
62 #include "names.h"
63 #include "force.h"
64 #include "disre.h"
65 #include "orires.h"
66 #include "pme.h"
67 #include "mdatoms.h"
68 #include "repl_ex.h"
69 #include "qmmm.h"
70 #include "domdec.h"
71 #include "domdec_network.h"
72 #include "partdec.h"
73 #include "topsort.h"
74 #include "coulomb.h"
75 #include "constr.h"
76 #include "shellfc.h"
77 #include "compute_io.h"
78 #include "mvdata.h"
79 #include "checkpoint.h"
80 #include "mtop_util.h"
81 #include "sighandler.h"
82 #include "txtdump.h"
83 #include "string2.h"
84 #include "pme_loadbal.h"
85 #include "bondf.h"
86 #include "membed.h"
87 #include "types/nlistheuristics.h"
88 #include "types/iteratedconstraints.h"
89 #include "nbnxn_cuda_data_mgmt.h"
90
91 #include "gromacs/utility/gmxmpi.h"
92
93 #ifdef GMX_FAHCORE
94 #include "corewrap.h"
95 #endif
96
97 static void reset_all_counters(FILE *fplog, t_commrec *cr,
98                                gmx_large_int_t step,
99                                gmx_large_int_t *step_rel, t_inputrec *ir,
100                                gmx_wallcycle_t wcycle, t_nrnb *nrnb,
101                                gmx_runtime_t *runtime,
102                                nbnxn_cuda_ptr_t cu_nbv)
103 {
104     char sbuf[STEPSTRSIZE];
105
106     /* Reset all the counters related to performance over the run */
107     md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
108                   gmx_step_str(step, sbuf));
109
110     if (cu_nbv)
111     {
112         nbnxn_cuda_reset_timings(cu_nbv);
113     }
114
115     wallcycle_stop(wcycle, ewcRUN);
116     wallcycle_reset_all(wcycle);
117     if (DOMAINDECOMP(cr))
118     {
119         reset_dd_statistics_counters(cr->dd);
120     }
121     init_nrnb(nrnb);
122     ir->init_step += *step_rel;
123     ir->nsteps    -= *step_rel;
124     *step_rel      = 0;
125     wallcycle_start(wcycle, ewcRUN);
126     runtime_start(runtime);
127     print_date_and_time(fplog, cr->nodeid, "Restarted time", runtime);
128 }
129
130 double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
131              const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
132              int nstglobalcomm,
133              gmx_vsite_t *vsite, gmx_constr_t constr,
134              int stepout, t_inputrec *ir,
135              gmx_mtop_t *top_global,
136              t_fcdata *fcd,
137              t_state *state_global,
138              t_mdatoms *mdatoms,
139              t_nrnb *nrnb, gmx_wallcycle_t wcycle,
140              gmx_edsam_t ed, t_forcerec *fr,
141              int repl_ex_nst, int repl_ex_nex, int repl_ex_seed, gmx_membed_t membed,
142              real cpt_period, real max_hours,
143              const char gmx_unused *deviceOptions,
144              unsigned long Flags,
145              gmx_runtime_t *runtime)
146 {
147     gmx_mdoutf_t   *outf;
148     gmx_large_int_t step, step_rel;
149     double          run_time;
150     double          t, t0, lam0[efptNR];
151     gmx_bool        bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
152     gmx_bool        bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
153                     bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
154                     bBornRadii, bStartingFromCpt;
155     gmx_bool          bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
156     gmx_bool          do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
157                       bForceUpdate = FALSE, bCPT;
158     gmx_bool          bMasterState;
159     int               force_flags, cglo_flags;
160     tensor            force_vir, shake_vir, total_vir, tmp_vir, pres;
161     int               i, m;
162     t_trxstatus      *status;
163     rvec              mu_tot;
164     t_vcm            *vcm;
165     t_state          *bufstate = NULL;
166     matrix           *scale_tot, pcoupl_mu, M, ebox;
167     gmx_nlheur_t      nlh;
168     t_trxframe        rerun_fr;
169     gmx_repl_ex_t     repl_ex = NULL;
170     int               nchkpt  = 1;
171     gmx_localtop_t   *top;
172     t_mdebin         *mdebin = NULL;
173     t_state          *state    = NULL;
174     rvec             *f_global = NULL;
175     gmx_enerdata_t   *enerd;
176     rvec             *f = NULL;
177     gmx_global_stat_t gstat;
178     gmx_update_t      upd   = NULL;
179     t_graph          *graph = NULL;
180     globsig_t         gs;
181     gmx_rng_t         mcrng = NULL;
182     gmx_groups_t     *groups;
183     gmx_ekindata_t   *ekind, *ekind_save;
184     gmx_shellfc_t     shellfc;
185     int               count, nconverged = 0;
186     real              timestep = 0;
187     double            tcount   = 0;
188     gmx_bool          bConverged = TRUE, bOK, bSumEkinhOld, bExchanged;
189     gmx_bool          bAppend;
190     gmx_bool          bResetCountersHalfMaxH = FALSE;
191     gmx_bool          bVV, bIterativeCase, bFirstIterate, bTemp, bPres, bTrotter;
192     gmx_bool          bUpdateDoLR;
193     real              dvdl_constr;
194     int               a0, a1;
195     rvec             *cbuf = NULL;
196     matrix            lastbox;
197     real              veta_save, scalevir, tracevir;
198     real              vetanew = 0;
199     int               lamnew  = 0;
200     /* for FEP */
201     int               nstfep;
202     double            cycles;
203     real              saved_conserved_quantity = 0;
204     real              last_ekin                = 0;
205     int               iter_i;
206     t_extmass         MassQ;
207     int             **trotter_seq;
208     char              sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
209     int               handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
210     gmx_iterate_t     iterate;
211     gmx_large_int_t   multisim_nsteps = -1;                        /* number of steps to do  before first multisim
212                                                                       simulation stops. If equal to zero, don't
213                                                                       communicate any more between multisims.*/
214     /* PME load balancing data for GPU kernels */
215     pme_load_balancing_t pme_loadbal = NULL;
216     double               cycles_pmes;
217     gmx_bool             bPMETuneTry = FALSE, bPMETuneRunning = FALSE;
218
219 #ifdef GMX_FAHCORE
220     /* Temporary addition for FAHCORE checkpointing */
221     int chkpt_ret;
222 #endif
223
224     /* Check for special mdrun options */
225     bRerunMD = (Flags & MD_RERUN);
226     bAppend  = (Flags & MD_APPENDFILES);
227     if (Flags & MD_RESETCOUNTERSHALFWAY)
228     {
229         if (ir->nsteps > 0)
230         {
231             /* Signal to reset the counters half the simulation steps. */
232             wcycle_set_reset_counters(wcycle, ir->nsteps/2);
233         }
234         /* Signal to reset the counters halfway the simulation time. */
235         bResetCountersHalfMaxH = (max_hours > 0);
236     }
237
238     /* md-vv uses averaged full step velocities for T-control
239        md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
240        md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
241     bVV = EI_VV(ir->eI);
242     if (bVV) /* to store the initial velocities while computing virial */
243     {
244         snew(cbuf, top_global->natoms);
245     }
246     /* all the iteratative cases - only if there are constraints */
247     bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
248     gmx_iterate_init(&iterate, FALSE); /* The default value of iterate->bIterationActive is set to
249                                           false in this step.  The correct value, true or false,
250                                           is set at each step, as it depends on the frequency of temperature
251                                           and pressure control.*/
252     bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
253
254     if (bRerunMD)
255     {
256         /* Since we don't know if the frames read are related in any way,
257          * rebuild the neighborlist at every step.
258          */
259         ir->nstlist       = 1;
260         ir->nstcalcenergy = 1;
261         nstglobalcomm     = 1;
262     }
263
264     check_ir_old_tpx_versions(cr, fplog, ir, top_global);
265
266     nstglobalcomm   = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
267     bGStatEveryStep = (nstglobalcomm == 1);
268
269     if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
270     {
271         fprintf(fplog,
272                 "To reduce the energy communication with nstlist = -1\n"
273                 "the neighbor list validity should not be checked at every step,\n"
274                 "this means that exact integration is not guaranteed.\n"
275                 "The neighbor list validity is checked after:\n"
276                 "  <n.list life time> - 2*std.dev.(n.list life time)  steps.\n"
277                 "In most cases this will result in exact integration.\n"
278                 "This reduces the energy communication by a factor of 2 to 3.\n"
279                 "If you want less energy communication, set nstlist > 3.\n\n");
280     }
281
282     if (bRerunMD)
283     {
284         ir->nstxtcout = 0;
285     }
286     groups = &top_global->groups;
287
288     /* Initial values */
289     init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
290             &(state_global->fep_state), lam0,
291             nrnb, top_global, &upd,
292             nfile, fnm, &outf, &mdebin,
293             force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags);
294
295     clear_mat(total_vir);
296     clear_mat(pres);
297     /* Energy terms and groups */
298     snew(enerd, 1);
299     init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
300                   enerd);
301     if (DOMAINDECOMP(cr))
302     {
303         f = NULL;
304     }
305     else
306     {
307         snew(f, top_global->natoms);
308     }
309
310     /* Kinetic energy data */
311     snew(ekind, 1);
312     init_ekindata(fplog, top_global, &(ir->opts), ekind);
313     /* needed for iteration of constraints */
314     snew(ekind_save, 1);
315     init_ekindata(fplog, top_global, &(ir->opts), ekind_save);
316     /* Copy the cos acceleration to the groups struct */
317     ekind->cosacc.cos_accel = ir->cos_accel;
318
319     gstat = global_stat_init(ir);
320     debug_gmx();
321
322     /* Check for polarizable models and flexible constraints */
323     shellfc = init_shell_flexcon(fplog,
324                                  top_global, n_flexible_constraints(constr),
325                                  (ir->bContinuation ||
326                                   (DOMAINDECOMP(cr) && !MASTER(cr))) ?
327                                  NULL : state_global->x);
328
329     if (DEFORM(*ir))
330     {
331 #ifdef GMX_THREAD_MPI
332         tMPI_Thread_mutex_lock(&deform_init_box_mutex);
333 #endif
334         set_deform_reference_box(upd,
335                                  deform_init_init_step_tpx,
336                                  deform_init_box_tpx);
337 #ifdef GMX_THREAD_MPI
338         tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
339 #endif
340     }
341
342     {
343         double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
344         if ((io > 2000) && MASTER(cr))
345         {
346             fprintf(stderr,
347                     "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
348                     io);
349         }
350     }
351
352     if (DOMAINDECOMP(cr))
353     {
354         top = dd_init_local_top(top_global);
355
356         snew(state, 1);
357         dd_init_local_state(cr->dd, state_global, state);
358
359         if (DDMASTER(cr->dd) && ir->nstfout)
360         {
361             snew(f_global, state_global->natoms);
362         }
363     }
364     else
365     {
366         if (PAR(cr))
367         {
368             /* Initialize the particle decomposition and split the topology */
369             top = split_system(fplog, top_global, ir, cr);
370
371             pd_cg_range(cr, &fr->cg0, &fr->hcg);
372             pd_at_range(cr, &a0, &a1);
373         }
374         else
375         {
376             top = gmx_mtop_generate_local_top(top_global, ir);
377
378             a0 = 0;
379             a1 = top_global->natoms;
380         }
381
382         forcerec_set_excl_load(fr, top, cr);
383
384         state    = partdec_init_local_state(cr, state_global);
385         f_global = f;
386
387         atoms2md(top_global, ir, 0, NULL, a0, a1-a0, mdatoms);
388
389         if (vsite)
390         {
391             set_vsite_top(vsite, top, mdatoms, cr);
392         }
393
394         if (ir->ePBC != epbcNONE && !fr->bMolPBC)
395         {
396             graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
397         }
398
399         if (shellfc)
400         {
401             make_local_shells(cr, mdatoms, shellfc);
402         }
403
404         setup_bonded_threading(fr, &top->idef);
405
406         if (ir->pull && PAR(cr))
407         {
408             dd_make_local_pull_groups(NULL, ir->pull, mdatoms);
409         }
410     }
411
412     if (DOMAINDECOMP(cr))
413     {
414         /* Distribute the charge groups over the nodes from the master node */
415         dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
416                             state_global, top_global, ir,
417                             state, &f, mdatoms, top, fr,
418                             vsite, shellfc, constr,
419                             nrnb, wcycle, FALSE);
420
421     }
422
423     update_mdatoms(mdatoms, state->lambda[efptMASS]);
424
425     if (opt2bSet("-cpi", nfile, fnm))
426     {
427         bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
428     }
429     else
430     {
431         bStateFromCP = FALSE;
432     }
433
434     if (ir->bExpanded)
435     {
436         init_expanded_ensemble(bStateFromCP,ir,&mcrng,&state->dfhist);
437     }
438
439     if (MASTER(cr))
440     {
441         if (bStateFromCP)
442         {
443             /* Update mdebin with energy history if appending to output files */
444             if (Flags & MD_APPENDFILES)
445             {
446                 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
447             }
448             else
449             {
450                 /* We might have read an energy history from checkpoint,
451                  * free the allocated memory and reset the counts.
452                  */
453                 done_energyhistory(&state_global->enerhist);
454                 init_energyhistory(&state_global->enerhist);
455             }
456         }
457         /* Set the initial energy history in state by updating once */
458         update_energyhistory(&state_global->enerhist, mdebin);
459     }
460
461     if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG))
462     {
463         /* Set the random state if we read a checkpoint file */
464         set_stochd_state(upd, state);
465     }
466
467     if (state->flags & (1<<estMC_RNG))
468     {
469         set_mc_state(mcrng, state);
470     }
471
472     /* Initialize constraints */
473     if (constr)
474     {
475         if (!DOMAINDECOMP(cr))
476         {
477             set_constraints(constr, top, ir, mdatoms, cr);
478         }
479     }
480
481     if (repl_ex_nst > 0)
482     {
483         /* We need to be sure replica exchange can only occur
484          * when the energies are current */
485         check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
486                         "repl_ex_nst", &repl_ex_nst);
487         /* This check needs to happen before inter-simulation
488          * signals are initialized, too */
489     }
490     if (repl_ex_nst > 0 && MASTER(cr))
491     {
492         repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
493                                         repl_ex_nst, repl_ex_nex, repl_ex_seed);
494     }
495
496     /* PME tuning is only supported with GPUs or PME nodes and not with rerun.
497      * With perturbed charges with soft-core we should not change the cut-off.
498      */
499     if ((Flags & MD_TUNEPME) &&
500         EEL_PME(fr->eeltype) &&
501         ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
502         !(ir->efep != efepNO && mdatoms->nChargePerturbed > 0 && ir->fepvals->bScCoul) &&
503         !bRerunMD)
504     {
505         pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
506         cycles_pmes = 0;
507         if (cr->duty & DUTY_PME)
508         {
509             /* Start tuning right away, as we can't measure the load */
510             bPMETuneRunning = TRUE;
511         }
512         else
513         {
514             /* Separate PME nodes, we can measure the PP/PME load balance */
515             bPMETuneTry = TRUE;
516         }
517     }
518
519     if (!ir->bContinuation && !bRerunMD)
520     {
521         if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
522         {
523             /* Set the velocities of frozen particles to zero */
524             for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
525             {
526                 for (m = 0; m < DIM; m++)
527                 {
528                     if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
529                     {
530                         state->v[i][m] = 0;
531                     }
532                 }
533             }
534         }
535
536         if (constr)
537         {
538             /* Constrain the initial coordinates and velocities */
539             do_constrain_first(fplog, constr, ir, mdatoms, state,
540                                cr, nrnb, fr, top);
541         }
542         if (vsite)
543         {
544             /* Construct the virtual sites for the initial configuration */
545             construct_vsites(vsite, state->x, ir->delta_t, NULL,
546                              top->idef.iparams, top->idef.il,
547                              fr->ePBC, fr->bMolPBC, graph, cr, state->box);
548         }
549     }
550
551     debug_gmx();
552
553     /* set free energy calculation frequency as the minimum
554        greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
555     nstfep = ir->fepvals->nstdhdl;
556     if (ir->bExpanded)
557     {
558         nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
559     }
560     if (repl_ex_nst > 0)
561     {
562         nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
563     }
564
565     /* I'm assuming we need global communication the first time! MRS */
566     cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
567                   | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
568                   | (bVV ? CGLO_PRESSURE : 0)
569                   | (bVV ? CGLO_CONSTRAINT : 0)
570                   | (bRerunMD ? CGLO_RERUNMD : 0)
571                   | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
572
573     bSumEkinhOld = FALSE;
574     compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
575                     NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
576                     constr, NULL, FALSE, state->box,
577                     top_global, &bSumEkinhOld, cglo_flags);
578     if (ir->eI == eiVVAK)
579     {
580         /* a second call to get the half step temperature initialized as well */
581         /* we do the same call as above, but turn the pressure off -- internally to
582            compute_globals, this is recognized as a velocity verlet half-step
583            kinetic energy calculation.  This minimized excess variables, but
584            perhaps loses some logic?*/
585
586         compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
587                         NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
588                         constr, NULL, FALSE, state->box,
589                         top_global, &bSumEkinhOld,
590                         cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
591     }
592
593     /* Calculate the initial half step temperature, and save the ekinh_old */
594     if (!(Flags & MD_STARTFROMCPT))
595     {
596         for (i = 0; (i < ir->opts.ngtc); i++)
597         {
598             copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
599         }
600     }
601     if (ir->eI != eiVV)
602     {
603         enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
604                                      and there is no previous step */
605     }
606
607     /* if using an iterative algorithm, we need to create a working directory for the state. */
608     if (bIterativeCase)
609     {
610         bufstate = init_bufstate(state);
611     }
612
613     /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
614        temperature control */
615     trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
616
617     if (MASTER(cr))
618     {
619         if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
620         {
621             fprintf(fplog,
622                     "RMS relative constraint deviation after constraining: %.2e\n",
623                     constr_rmsd(constr, FALSE));
624         }
625         if (EI_STATE_VELOCITY(ir->eI))
626         {
627             fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
628         }
629         if (bRerunMD)
630         {
631             fprintf(stderr, "starting md rerun '%s', reading coordinates from"
632                     " input trajectory '%s'\n\n",
633                     *(top_global->name), opt2fn("-rerun", nfile, fnm));
634             if (bVerbose)
635             {
636                 fprintf(stderr, "Calculated time to finish depends on nsteps from "
637                         "run input file,\nwhich may not correspond to the time "
638                         "needed to process input trajectory.\n\n");
639             }
640         }
641         else
642         {
643             char tbuf[20];
644             fprintf(stderr, "starting mdrun '%s'\n",
645                     *(top_global->name));
646             if (ir->nsteps >= 0)
647             {
648                 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
649             }
650             else
651             {
652                 sprintf(tbuf, "%s", "infinite");
653             }
654             if (ir->init_step > 0)
655             {
656                 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
657                         gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
658                         gmx_step_str(ir->init_step, sbuf2),
659                         ir->init_step*ir->delta_t);
660             }
661             else
662             {
663                 fprintf(stderr, "%s steps, %s ps.\n",
664                         gmx_step_str(ir->nsteps, sbuf), tbuf);
665             }
666         }
667         fprintf(fplog, "\n");
668     }
669
670     /* Set and write start time */
671     runtime_start(runtime);
672     print_date_and_time(fplog, cr->nodeid, "Started mdrun", runtime);
673     wallcycle_start(wcycle, ewcRUN);
674     if (fplog)
675     {
676         fprintf(fplog, "\n");
677     }
678
679     /* safest point to do file checkpointing is here.  More general point would be immediately before integrator call */
680 #ifdef GMX_FAHCORE
681     chkpt_ret = fcCheckPointParallel( cr->nodeid,
682                                       NULL, 0);
683     if (chkpt_ret == 0)
684     {
685         gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
686     }
687 #endif
688
689     debug_gmx();
690     /***********************************************************
691      *
692      *             Loop over MD steps
693      *
694      ************************************************************/
695
696     /* if rerunMD then read coordinates and velocities from input trajectory */
697     if (bRerunMD)
698     {
699         if (getenv("GMX_FORCE_UPDATE"))
700         {
701             bForceUpdate = TRUE;
702         }
703
704         rerun_fr.natoms = 0;
705         if (MASTER(cr))
706         {
707             bNotLastFrame = read_first_frame(oenv, &status,
708                                              opt2fn("-rerun", nfile, fnm),
709                                              &rerun_fr, TRX_NEED_X | TRX_READ_V);
710             if (rerun_fr.natoms != top_global->natoms)
711             {
712                 gmx_fatal(FARGS,
713                           "Number of atoms in trajectory (%d) does not match the "
714                           "run input file (%d)\n",
715                           rerun_fr.natoms, top_global->natoms);
716             }
717             if (ir->ePBC != epbcNONE)
718             {
719                 if (!rerun_fr.bBox)
720                 {
721                     gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f does not contain a box, while pbc is used", rerun_fr.step, rerun_fr.time);
722                 }
723                 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
724                 {
725                     gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
726                 }
727             }
728         }
729
730         if (PAR(cr))
731         {
732             rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
733         }
734
735         if (ir->ePBC != epbcNONE)
736         {
737             /* Set the shift vectors.
738              * Necessary here when have a static box different from the tpr box.
739              */
740             calc_shifts(rerun_fr.box, fr->shift_vec);
741         }
742     }
743
744     /* loop over MD steps or if rerunMD to end of input trajectory */
745     bFirstStep = TRUE;
746     /* Skip the first Nose-Hoover integration when we get the state from tpx */
747     bStateFromTPX    = !bStateFromCP;
748     bInitStep        = bFirstStep && (bStateFromTPX || bVV);
749     bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
750     bLastStep        = FALSE;
751     bSumEkinhOld     = FALSE;
752     bExchanged       = FALSE;
753
754     init_global_signals(&gs, cr, ir, repl_ex_nst);
755
756     step     = ir->init_step;
757     step_rel = 0;
758
759     if (ir->nstlist == -1)
760     {
761         init_nlistheuristics(&nlh, bGStatEveryStep, step);
762     }
763
764     if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
765     {
766         /* check how many steps are left in other sims */
767         multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
768     }
769
770
771     /* and stop now if we should */
772     bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
773                  ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
774     while (!bLastStep || (bRerunMD && bNotLastFrame))
775     {
776
777         wallcycle_start(wcycle, ewcSTEP);
778
779         if (bRerunMD)
780         {
781             if (rerun_fr.bStep)
782             {
783                 step     = rerun_fr.step;
784                 step_rel = step - ir->init_step;
785             }
786             if (rerun_fr.bTime)
787             {
788                 t = rerun_fr.time;
789             }
790             else
791             {
792                 t = step;
793             }
794         }
795         else
796         {
797             bLastStep = (step_rel == ir->nsteps);
798             t         = t0 + step*ir->delta_t;
799         }
800
801         if (ir->efep != efepNO || ir->bSimTemp)
802         {
803             /* find and set the current lambdas.  If rerunning, we either read in a state, or a lambda value,
804                requiring different logic. */
805
806             set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
807             bDoDHDL      = do_per_step(step, ir->fepvals->nstdhdl);
808             bDoFEP       = (do_per_step(step, nstfep) && (ir->efep != efepNO));
809             bDoExpanded  = (do_per_step(step, ir->expandedvals->nstexpanded)
810                             && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
811         }
812
813         if (bSimAnn)
814         {
815             update_annealing_target_temp(&(ir->opts), t);
816         }
817
818         if (bRerunMD)
819         {
820             if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
821             {
822                 for (i = 0; i < state_global->natoms; i++)
823                 {
824                     copy_rvec(rerun_fr.x[i], state_global->x[i]);
825                 }
826                 if (rerun_fr.bV)
827                 {
828                     for (i = 0; i < state_global->natoms; i++)
829                     {
830                         copy_rvec(rerun_fr.v[i], state_global->v[i]);
831                     }
832                 }
833                 else
834                 {
835                     for (i = 0; i < state_global->natoms; i++)
836                     {
837                         clear_rvec(state_global->v[i]);
838                     }
839                     if (bRerunWarnNoV)
840                     {
841                         fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
842                                 "         Ekin, temperature and pressure are incorrect,\n"
843                                 "         the virial will be incorrect when constraints are present.\n"
844                                 "\n");
845                         bRerunWarnNoV = FALSE;
846                     }
847                 }
848             }
849             copy_mat(rerun_fr.box, state_global->box);
850             copy_mat(state_global->box, state->box);
851
852             if (vsite && (Flags & MD_RERUN_VSITE))
853             {
854                 if (DOMAINDECOMP(cr))
855                 {
856                     gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
857                 }
858                 if (graph)
859                 {
860                     /* Following is necessary because the graph may get out of sync
861                      * with the coordinates if we only have every N'th coordinate set
862                      */
863                     mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
864                     shift_self(graph, state->box, state->x);
865                 }
866                 construct_vsites(vsite, state->x, ir->delta_t, state->v,
867                                  top->idef.iparams, top->idef.il,
868                                  fr->ePBC, fr->bMolPBC, graph, cr, state->box);
869                 if (graph)
870                 {
871                     unshift_self(graph, state->box, state->x);
872                 }
873             }
874         }
875
876         /* Stop Center of Mass motion */
877         bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
878
879         if (bRerunMD)
880         {
881             /* for rerun MD always do Neighbour Searching */
882             bNS      = (bFirstStep || ir->nstlist != 0);
883             bNStList = bNS;
884         }
885         else
886         {
887             /* Determine whether or not to do Neighbour Searching and LR */
888             bNStList = (ir->nstlist > 0  && step % ir->nstlist == 0);
889
890             bNS = (bFirstStep || bExchanged || bNStList || bDoFEP ||
891                    (ir->nstlist == -1 && nlh.nabnsb > 0));
892
893             if (bNS && ir->nstlist == -1)
894             {
895                 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bDoFEP, step);
896             }
897         }
898
899         /* check whether we should stop because another simulation has
900            stopped. */
901         if (MULTISIM(cr))
902         {
903             if ( (multisim_nsteps >= 0) &&  (step_rel >= multisim_nsteps)  &&
904                  (multisim_nsteps != ir->nsteps) )
905             {
906                 if (bNS)
907                 {
908                     if (MASTER(cr))
909                     {
910                         fprintf(stderr,
911                                 "Stopping simulation %d because another one has finished\n",
912                                 cr->ms->sim);
913                     }
914                     bLastStep         = TRUE;
915                     gs.sig[eglsCHKPT] = 1;
916                 }
917             }
918         }
919
920         /* < 0 means stop at next step, > 0 means stop at next NS step */
921         if ( (gs.set[eglsSTOPCOND] < 0) ||
922              ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
923         {
924             bLastStep = TRUE;
925         }
926
927         /* Determine whether or not to update the Born radii if doing GB */
928         bBornRadii = bFirstStep;
929         if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
930         {
931             bBornRadii = TRUE;
932         }
933
934         do_log     = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
935         do_verbose = bVerbose &&
936             (step % stepout == 0 || bFirstStep || bLastStep);
937
938         if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
939         {
940             if (bRerunMD)
941             {
942                 bMasterState = TRUE;
943             }
944             else
945             {
946                 bMasterState = FALSE;
947                 /* Correct the new box if it is too skewed */
948                 if (DYNAMIC_BOX(*ir))
949                 {
950                     if (correct_box(fplog, step, state->box, graph))
951                     {
952                         bMasterState = TRUE;
953                     }
954                 }
955                 if (DOMAINDECOMP(cr) && bMasterState)
956                 {
957                     dd_collect_state(cr->dd, state, state_global);
958                 }
959             }
960
961             if (DOMAINDECOMP(cr))
962             {
963                 /* Repartition the domain decomposition */
964                 wallcycle_start(wcycle, ewcDOMDEC);
965                 dd_partition_system(fplog, step, cr,
966                                     bMasterState, nstglobalcomm,
967                                     state_global, top_global, ir,
968                                     state, &f, mdatoms, top, fr,
969                                     vsite, shellfc, constr,
970                                     nrnb, wcycle,
971                                     do_verbose && !bPMETuneRunning);
972                 wallcycle_stop(wcycle, ewcDOMDEC);
973                 /* If using an iterative integrator, reallocate space to match the decomposition */
974             }
975         }
976
977         if (MASTER(cr) && do_log)
978         {
979             print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
980         }
981
982         if (ir->efep != efepNO)
983         {
984             update_mdatoms(mdatoms, state->lambda[efptMASS]);
985         }
986
987         if ((bRerunMD && rerun_fr.bV) || bExchanged)
988         {
989
990             /* We need the kinetic energy at minus the half step for determining
991              * the full step kinetic energy and possibly for T-coupling.*/
992             /* This may not be quite working correctly yet . . . . */
993             compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
994                             wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
995                             constr, NULL, FALSE, state->box,
996                             top_global, &bSumEkinhOld,
997                             CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
998         }
999         clear_mat(force_vir);
1000
1001         /* We write a checkpoint at this MD step when:
1002          * either at an NS step when we signalled through gs,
1003          * or at the last step (but not when we do not want confout),
1004          * but never at the first step or with rerun.
1005          */
1006         bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1007                  (bLastStep && (Flags & MD_CONFOUT))) &&
1008                 step > ir->init_step && !bRerunMD);
1009         if (bCPT)
1010         {
1011             gs.set[eglsCHKPT] = 0;
1012         }
1013
1014         /* Determine the energy and pressure:
1015          * at nstcalcenergy steps and at energy output steps (set below).
1016          */
1017         if (EI_VV(ir->eI) && (!bInitStep))
1018         {
1019             /* for vv, the first half of the integration actually corresponds
1020                to the previous step.  bCalcEner is only required to be evaluated on the 'next' step,
1021                but the virial needs to be calculated on both the current step and the 'next' step. Future
1022                reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1023
1024             bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1025             bCalcVir  = bCalcEner ||
1026                 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1027         }
1028         else
1029         {
1030             bCalcEner = do_per_step(step, ir->nstcalcenergy);
1031             bCalcVir  = bCalcEner ||
1032                 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1033         }
1034
1035         /* Do we need global communication ? */
1036         bGStat = (bCalcVir || bCalcEner || bStopCM ||
1037                   do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1038                   (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1039
1040         do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1041
1042         if (do_ene || do_log)
1043         {
1044             bCalcVir  = TRUE;
1045             bCalcEner = TRUE;
1046             bGStat    = TRUE;
1047         }
1048
1049         /* these CGLO_ options remain the same throughout the iteration */
1050         cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1051                       (bGStat ? CGLO_GSTAT : 0)
1052                       );
1053
1054         force_flags = (GMX_FORCE_STATECHANGED |
1055                        ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1056                        GMX_FORCE_ALLFORCES |
1057                        GMX_FORCE_SEPLRF |
1058                        (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1059                        (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1060                        (bDoFEP ? GMX_FORCE_DHDL : 0)
1061                        );
1062
1063         if (fr->bTwinRange)
1064         {
1065             if (do_per_step(step, ir->nstcalclr))
1066             {
1067                 force_flags |= GMX_FORCE_DO_LR;
1068             }
1069         }
1070
1071         if (shellfc)
1072         {
1073             /* Now is the time to relax the shells */
1074             count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1075                                         ir, bNS, force_flags,
1076                                         top,
1077                                         constr, enerd, fcd,
1078                                         state, f, force_vir, mdatoms,
1079                                         nrnb, wcycle, graph, groups,
1080                                         shellfc, fr, bBornRadii, t, mu_tot,
1081                                         &bConverged, vsite,
1082                                         outf->fp_field);
1083             tcount += count;
1084
1085             if (bConverged)
1086             {
1087                 nconverged++;
1088             }
1089         }
1090         else
1091         {
1092             /* The coordinates (x) are shifted (to get whole molecules)
1093              * in do_force.
1094              * This is parallellized as well, and does communication too.
1095              * Check comments in sim_util.c
1096              */
1097             do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1098                      state->box, state->x, &state->hist,
1099                      f, force_vir, mdatoms, enerd, fcd,
1100                      state->lambda, graph,
1101                      fr, vsite, mu_tot, t, outf->fp_field, ed, bBornRadii,
1102                      (bNS ? GMX_FORCE_NS : 0) | force_flags);
1103         }
1104
1105         if (bVV && !bStartingFromCpt && !bRerunMD)
1106         /*  ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1107         {
1108             if (ir->eI == eiVV && bInitStep)
1109             {
1110                 /* if using velocity verlet with full time step Ekin,
1111                  * take the first half step only to compute the
1112                  * virial for the first step. From there,
1113                  * revert back to the initial coordinates
1114                  * so that the input is actually the initial step.
1115                  */
1116                 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1117             }
1118             else
1119             {
1120                 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1121                 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1122             }
1123
1124             /* If we are using twin-range interactions where the long-range component
1125              * is only evaluated every nstcalclr>1 steps, we should do a special update
1126              * step to combine the long-range forces on these steps.
1127              * For nstcalclr=1 this is not done, since the forces would have been added
1128              * directly to the short-range forces already.
1129              */
1130             bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1131
1132             update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1133                           f, bUpdateDoLR, fr->f_twin, fcd,
1134                           ekind, M, upd, bInitStep, etrtVELOCITY1,
1135                           cr, nrnb, constr, &top->idef);
1136
1137             if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1138             {
1139                 gmx_iterate_init(&iterate, TRUE);
1140             }
1141             /* for iterations, we save these vectors, as we will be self-consistently iterating
1142                the calculations */
1143
1144             /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1145
1146             /* save the state */
1147             if (iterate.bIterationActive)
1148             {
1149                 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1150             }
1151
1152             bFirstIterate = TRUE;
1153             while (bFirstIterate || iterate.bIterationActive)
1154             {
1155                 if (iterate.bIterationActive)
1156                 {
1157                     copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1158                     if (bFirstIterate && bTrotter)
1159                     {
1160                         /* The first time through, we need a decent first estimate
1161                            of veta(t+dt) to compute the constraints.  Do
1162                            this by computing the box volume part of the
1163                            trotter integration at this time. Nothing else
1164                            should be changed by this routine here.  If
1165                            !(first time), we start with the previous value
1166                            of veta.  */
1167
1168                         veta_save = state->veta;
1169                         trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1170                         vetanew     = state->veta;
1171                         state->veta = veta_save;
1172                     }
1173                 }
1174
1175                 bOK = TRUE;
1176                 if (!bRerunMD || rerun_fr.bV || bForceUpdate)     /* Why is rerun_fr.bV here?  Unclear. */
1177                 {
1178                     update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1179                                        state, fr->bMolPBC, graph, f,
1180                                        &top->idef, shake_vir,
1181                                        cr, nrnb, wcycle, upd, constr,
1182                                        TRUE, bCalcVir, vetanew);
1183
1184                     if (!bOK)
1185                     {
1186                         gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1187                     }
1188
1189                 }
1190                 else if (graph)
1191                 {
1192                     /* Need to unshift here if a do_force has been
1193                        called in the previous step */
1194                     unshift_self(graph, state->box, state->x);
1195                 }
1196
1197                 /* if VV, compute the pressure and constraints */
1198                 /* For VV2, we strictly only need this if using pressure
1199                  * control, but we really would like to have accurate pressures
1200                  * printed out.
1201                  * Think about ways around this in the future?
1202                  * For now, keep this choice in comments.
1203                  */
1204                 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1205                 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1206                 bPres = TRUE;
1207                 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1208                 if (bCalcEner && ir->eI == eiVVAK)  /*MRS:  7/9/2010 -- this still doesn't fix it?*/
1209                 {
1210                     bSumEkinhOld = TRUE;
1211                 }
1212                 /* for vv, the first half of the integration actually corresponds to the previous step.
1213                    So we need information from the last step in the first half of the integration */
1214                 if (bGStat || do_per_step(step-1, nstglobalcomm))
1215                 {
1216                     compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1217                                     wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1218                                     constr, NULL, FALSE, state->box,
1219                                     top_global, &bSumEkinhOld,
1220                                     cglo_flags
1221                                     | CGLO_ENERGY
1222                                     | (bTemp ? CGLO_TEMPERATURE : 0)
1223                                     | (bPres ? CGLO_PRESSURE : 0)
1224                                     | (bPres ? CGLO_CONSTRAINT : 0)
1225                                     | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1226                                     | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1227                                     | CGLO_SCALEEKIN
1228                                     );
1229                     /* explanation of above:
1230                        a) We compute Ekin at the full time step
1231                        if 1) we are using the AveVel Ekin, and it's not the
1232                        initial step, or 2) if we are using AveEkin, but need the full
1233                        time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1234                        b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1235                        EkinAveVel because it's needed for the pressure */
1236                 }
1237                 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1238                 if (!bInitStep)
1239                 {
1240                     if (bTrotter)
1241                     {
1242                         m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1243                         trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1244                     }
1245                     else
1246                     {
1247                         if (bExchanged)
1248                         {
1249
1250                             /* We need the kinetic energy at minus the half step for determining
1251                              * the full step kinetic energy and possibly for T-coupling.*/
1252                             /* This may not be quite working correctly yet . . . . */
1253                             compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1254                                             wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1255                                             constr, NULL, FALSE, state->box,
1256                                             top_global, &bSumEkinhOld,
1257                                             CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1258                         }
1259                     }
1260                 }
1261
1262                 if (iterate.bIterationActive &&
1263                     done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1264                                    state->veta, &vetanew))
1265                 {
1266                     break;
1267                 }
1268                 bFirstIterate = FALSE;
1269             }
1270
1271             if (bTrotter && !bInitStep)
1272             {
1273                 copy_mat(shake_vir, state->svir_prev);
1274                 copy_mat(force_vir, state->fvir_prev);
1275                 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1276                 {
1277                     /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1278                     enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1279                     enerd->term[F_EKIN] = trace(ekind->ekin);
1280                 }
1281             }
1282             /* if it's the initial step, we performed this first step just to get the constraint virial */
1283             if (bInitStep && ir->eI == eiVV)
1284             {
1285                 copy_rvecn(cbuf, state->v, 0, state->natoms);
1286             }
1287         }
1288
1289         /* MRS -- now done iterating -- compute the conserved quantity */
1290         if (bVV)
1291         {
1292             saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1293             if (ir->eI == eiVV)
1294             {
1295                 last_ekin = enerd->term[F_EKIN];
1296             }
1297             if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1298             {
1299                 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1300             }
1301             /* sum up the foreign energy and dhdl terms for vv.  currently done every step so that dhdl is correct in the .edr */
1302             if (!bRerunMD)
1303             {
1304                 sum_dhdl(enerd, state->lambda, ir->fepvals);
1305             }
1306         }
1307
1308         /* ########  END FIRST UPDATE STEP  ############## */
1309         /* ########  If doing VV, we now have v(dt) ###### */
1310         if (bDoExpanded)
1311         {
1312             /* perform extended ensemble sampling in lambda - we don't
1313                actually move to the new state before outputting
1314                statistics, but if performing simulated tempering, we
1315                do update the velocities and the tau_t. */
1316
1317             lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, mcrng, state->v, mdatoms);
1318             /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1319             copy_df_history(&state_global->dfhist,&state->dfhist);
1320         }
1321
1322         /* Now we have the energies and forces corresponding to the
1323          * coordinates at time t. We must output all of this before
1324          * the update.
1325          */
1326         do_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1327                               ir, state, state_global, top_global, fr, upd,
1328                               outf, mdebin, ekind, f, f_global,
1329                               wcycle, mcrng, &nchkpt,
1330                               bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1331                               bSumEkinhOld);
1332
1333         /* kludge -- virial is lost with restart for NPT control. Must restart */
1334         if (bStartingFromCpt && bVV)
1335         {
1336             copy_mat(state->svir_prev, shake_vir);
1337             copy_mat(state->fvir_prev, force_vir);
1338         }
1339
1340         /* Determine the wallclock run time up till now */
1341         run_time = gmx_gettime() - (double)runtime->real;
1342
1343         /* Check whether everything is still allright */
1344         if (((int)gmx_get_stop_condition() > handled_stop_condition)
1345 #ifdef GMX_THREAD_MPI
1346             && MASTER(cr)
1347 #endif
1348             )
1349         {
1350             /* this is just make gs.sig compatible with the hack
1351                of sending signals around by MPI_Reduce with together with
1352                other floats */
1353             if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1354             {
1355                 gs.sig[eglsSTOPCOND] = 1;
1356             }
1357             if (gmx_get_stop_condition() == gmx_stop_cond_next)
1358             {
1359                 gs.sig[eglsSTOPCOND] = -1;
1360             }
1361             /* < 0 means stop at next step, > 0 means stop at next NS step */
1362             if (fplog)
1363             {
1364                 fprintf(fplog,
1365                         "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1366                         gmx_get_signal_name(),
1367                         gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1368                 fflush(fplog);
1369             }
1370             fprintf(stderr,
1371                     "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1372                     gmx_get_signal_name(),
1373                     gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1374             fflush(stderr);
1375             handled_stop_condition = (int)gmx_get_stop_condition();
1376         }
1377         else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1378                  (max_hours > 0 && run_time > max_hours*60.0*60.0*0.99) &&
1379                  gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1380         {
1381             /* Signal to terminate the run */
1382             gs.sig[eglsSTOPCOND] = 1;
1383             if (fplog)
1384             {
1385                 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1386             }
1387             fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1388         }
1389
1390         if (bResetCountersHalfMaxH && MASTER(cr) &&
1391             run_time > max_hours*60.0*60.0*0.495)
1392         {
1393             gs.sig[eglsRESETCOUNTERS] = 1;
1394         }
1395
1396         if (ir->nstlist == -1 && !bRerunMD)
1397         {
1398             /* When bGStatEveryStep=FALSE, global_stat is only called
1399              * when we check the atom displacements, not at NS steps.
1400              * This means that also the bonded interaction count check is not
1401              * performed immediately after NS. Therefore a few MD steps could
1402              * be performed with missing interactions.
1403              * But wrong energies are never written to file,
1404              * since energies are only written after global_stat
1405              * has been called.
1406              */
1407             if (step >= nlh.step_nscheck)
1408             {
1409                 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1410                                                      nlh.scale_tot, state->x);
1411             }
1412             else
1413             {
1414                 /* This is not necessarily true,
1415                  * but step_nscheck is determined quite conservatively.
1416                  */
1417                 nlh.nabnsb = 0;
1418             }
1419         }
1420
1421         /* In parallel we only have to check for checkpointing in steps
1422          * where we do global communication,
1423          *  otherwise the other nodes don't know.
1424          */
1425         if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1426                            cpt_period >= 0 &&
1427                            (cpt_period == 0 ||
1428                             run_time >= nchkpt*cpt_period*60.0)) &&
1429             gs.set[eglsCHKPT] == 0)
1430         {
1431             gs.sig[eglsCHKPT] = 1;
1432         }
1433
1434         /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1435         if (EI_VV(ir->eI))
1436         {
1437             if (!bInitStep)
1438             {
1439                 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1440             }
1441             if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1442             {
1443                 gmx_bool bIfRandomize;
1444                 bIfRandomize = update_randomize_velocities(ir, step, mdatoms, state, upd, &top->idef, constr);
1445                 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1446                 if (constr && bIfRandomize)
1447                 {
1448                     update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1449                                        state, fr->bMolPBC, graph, f,
1450                                        &top->idef, tmp_vir,
1451                                        cr, nrnb, wcycle, upd, constr,
1452                                        TRUE, bCalcVir, vetanew);
1453                 }
1454             }
1455         }
1456
1457         if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1458         {
1459             gmx_iterate_init(&iterate, TRUE);
1460             /* for iterations, we save these vectors, as we will be redoing the calculations */
1461             copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1462         }
1463
1464         bFirstIterate = TRUE;
1465         while (bFirstIterate || iterate.bIterationActive)
1466         {
1467             /* We now restore these vectors to redo the calculation with improved extended variables */
1468             if (iterate.bIterationActive)
1469             {
1470                 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1471             }
1472
1473             /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1474                so scroll down for that logic */
1475
1476             /* #########   START SECOND UPDATE STEP ################# */
1477             /* Box is changed in update() when we do pressure coupling,
1478              * but we should still use the old box for energy corrections and when
1479              * writing it to the energy file, so it matches the trajectory files for
1480              * the same timestep above. Make a copy in a separate array.
1481              */
1482             copy_mat(state->box, lastbox);
1483
1484             bOK         = TRUE;
1485             dvdl_constr = 0;
1486
1487             if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1488             {
1489                 wallcycle_start(wcycle, ewcUPDATE);
1490                 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1491                 if (bTrotter)
1492                 {
1493                     if (iterate.bIterationActive)
1494                     {
1495                         if (bFirstIterate)
1496                         {
1497                             scalevir = 1;
1498                         }
1499                         else
1500                         {
1501                             /* we use a new value of scalevir to converge the iterations faster */
1502                             scalevir = tracevir/trace(shake_vir);
1503                         }
1504                         msmul(shake_vir, scalevir, shake_vir);
1505                         m_add(force_vir, shake_vir, total_vir);
1506                         clear_mat(shake_vir);
1507                     }
1508                     trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1509                     /* We can only do Berendsen coupling after we have summed
1510                      * the kinetic energy or virial. Since the happens
1511                      * in global_state after update, we should only do it at
1512                      * step % nstlist = 1 with bGStatEveryStep=FALSE.
1513                      */
1514                 }
1515                 else
1516                 {
1517                     update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1518                     update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1519                 }
1520
1521                 if (bVV)
1522                 {
1523                     bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1524
1525                     /* velocity half-step update */
1526                     update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1527                                   bUpdateDoLR, fr->f_twin, fcd,
1528                                   ekind, M, upd, FALSE, etrtVELOCITY2,
1529                                   cr, nrnb, constr, &top->idef);
1530                 }
1531
1532                 /* Above, initialize just copies ekinh into ekin,
1533                  * it doesn't copy position (for VV),
1534                  * and entire integrator for MD.
1535                  */
1536
1537                 if (ir->eI == eiVVAK)
1538                 {
1539                     copy_rvecn(state->x, cbuf, 0, state->natoms);
1540                 }
1541                 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1542
1543                 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1544                               bUpdateDoLR, fr->f_twin, fcd,
1545                               ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1546                 wallcycle_stop(wcycle, ewcUPDATE);
1547
1548                 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1549                                    fr->bMolPBC, graph, f,
1550                                    &top->idef, shake_vir,
1551                                    cr, nrnb, wcycle, upd, constr,
1552                                    FALSE, bCalcVir, state->veta);
1553
1554                 if (ir->eI == eiVVAK)
1555                 {
1556                     /* erase F_EKIN and F_TEMP here? */
1557                     /* just compute the kinetic energy at the half step to perform a trotter step */
1558                     compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1559                                     wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1560                                     constr, NULL, FALSE, lastbox,
1561                                     top_global, &bSumEkinhOld,
1562                                     cglo_flags | CGLO_TEMPERATURE
1563                                     );
1564                     wallcycle_start(wcycle, ewcUPDATE);
1565                     trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1566                     /* now we know the scaling, we can compute the positions again again */
1567                     copy_rvecn(cbuf, state->x, 0, state->natoms);
1568
1569                     bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1570
1571                     update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1572                                   bUpdateDoLR, fr->f_twin, fcd,
1573                                   ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1574                     wallcycle_stop(wcycle, ewcUPDATE);
1575
1576                     /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1577                     /* are the small terms in the shake_vir here due
1578                      * to numerical errors, or are they important
1579                      * physically? I'm thinking they are just errors, but not completely sure.
1580                      * For now, will call without actually constraining, constr=NULL*/
1581                     update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1582                                        state, fr->bMolPBC, graph, f,
1583                                        &top->idef, tmp_vir,
1584                                        cr, nrnb, wcycle, upd, NULL,
1585                                        FALSE, bCalcVir,
1586                                        state->veta);
1587                 }
1588                 if (!bOK)
1589                 {
1590                     gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1591                 }
1592
1593                 if (fr->bSepDVDL && fplog && do_log)
1594                 {
1595                     gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1596                 }
1597                 if (bVV)
1598                 {
1599                     /* this factor or 2 correction is necessary
1600                        because half of the constraint force is removed
1601                        in the vv step, so we have to double it.  See
1602                        the Redmine issue #1255.  It is not yet clear
1603                        if the factor of 2 is exact, or just a very
1604                        good approximation, and this will be
1605                        investigated.  The next step is to see if this
1606                        can be done adding a dhdl contribution from the
1607                        rattle step, but this is somewhat more
1608                        complicated with the current code. Will be
1609                        investigated, hopefully for 4.6.3. However,
1610                        this current solution is much better than
1611                        having it completely wrong.
1612                      */
1613                     enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1614                 }
1615                 else
1616                 {
1617                     enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1618                 }
1619             }
1620             else if (graph)
1621             {
1622                 /* Need to unshift here */
1623                 unshift_self(graph, state->box, state->x);
1624             }
1625
1626             if (vsite != NULL)
1627             {
1628                 wallcycle_start(wcycle, ewcVSITECONSTR);
1629                 if (graph != NULL)
1630                 {
1631                     shift_self(graph, state->box, state->x);
1632                 }
1633                 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1634                                  top->idef.iparams, top->idef.il,
1635                                  fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1636
1637                 if (graph != NULL)
1638                 {
1639                     unshift_self(graph, state->box, state->x);
1640                 }
1641                 wallcycle_stop(wcycle, ewcVSITECONSTR);
1642             }
1643
1644             /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints  ############ */
1645             /* With Leap-Frog we can skip compute_globals at
1646              * non-communication steps, but we need to calculate
1647              * the kinetic energy one step before communication.
1648              */
1649             if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1650             {
1651                 if (ir->nstlist == -1 && bFirstIterate)
1652                 {
1653                     gs.sig[eglsNABNSB] = nlh.nabnsb;
1654                 }
1655                 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1656                                 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1657                                 constr,
1658                                 bFirstIterate ? &gs : NULL,
1659                                 (step_rel % gs.nstms == 0) &&
1660                                 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1661                                 lastbox,
1662                                 top_global, &bSumEkinhOld,
1663                                 cglo_flags
1664                                 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1665                                 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1666                                 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1667                                 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1668                                 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1669                                 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1670                                 | CGLO_CONSTRAINT
1671                                 );
1672                 if (ir->nstlist == -1 && bFirstIterate)
1673                 {
1674                     nlh.nabnsb         = gs.set[eglsNABNSB];
1675                     gs.set[eglsNABNSB] = 0;
1676                 }
1677             }
1678             /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1679             /* #############  END CALC EKIN AND PRESSURE ################# */
1680
1681             /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1682                the virial that should probably be addressed eventually. state->veta has better properies,
1683                but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1684                generate the new shake_vir, but test the veta value for convergence.  This will take some thought. */
1685
1686             if (iterate.bIterationActive &&
1687                 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1688                                trace(shake_vir), &tracevir))
1689             {
1690                 break;
1691             }
1692             bFirstIterate = FALSE;
1693         }
1694
1695         if (!bVV || bRerunMD)
1696         {
1697             /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1698             sum_dhdl(enerd, state->lambda, ir->fepvals);
1699         }
1700         update_box(fplog, step, ir, mdatoms, state, f,
1701                    ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1702
1703         /* ################# END UPDATE STEP 2 ################# */
1704         /* #### We now have r(t+dt) and v(t+dt/2)  ############# */
1705
1706         /* The coordinates (x) were unshifted in update */
1707         if (!bGStat)
1708         {
1709             /* We will not sum ekinh_old,
1710              * so signal that we still have to do it.
1711              */
1712             bSumEkinhOld = TRUE;
1713         }
1714
1715         /* #########  BEGIN PREPARING EDR OUTPUT  ###########  */
1716
1717         /* use the directly determined last velocity, not actually the averaged half steps */
1718         if (bTrotter && ir->eI == eiVV)
1719         {
1720             enerd->term[F_EKIN] = last_ekin;
1721         }
1722         enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1723
1724         if (bVV)
1725         {
1726             enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1727         }
1728         else
1729         {
1730             enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1731         }
1732         /* #########  END PREPARING EDR OUTPUT  ###########  */
1733
1734         /* Time for performance */
1735         if (((step % stepout) == 0) || bLastStep)
1736         {
1737             runtime_upd_proc(runtime);
1738         }
1739
1740         /* Output stuff */
1741         if (MASTER(cr))
1742         {
1743             gmx_bool do_dr, do_or;
1744
1745             if (fplog && do_log && bDoExpanded)
1746             {
1747                 /* only needed if doing expanded ensemble */
1748                 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1749                                           &state_global->dfhist, state->fep_state, ir->nstlog, step);
1750             }
1751             if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1752             {
1753                 if (bCalcEner)
1754                 {
1755                     upd_mdebin(mdebin, bDoDHDL, TRUE,
1756                                t, mdatoms->tmass, enerd, state,
1757                                ir->fepvals, ir->expandedvals, lastbox,
1758                                shake_vir, force_vir, total_vir, pres,
1759                                ekind, mu_tot, constr);
1760                 }
1761                 else
1762                 {
1763                     upd_mdebin_step(mdebin);
1764                 }
1765
1766                 do_dr  = do_per_step(step, ir->nstdisreout);
1767                 do_or  = do_per_step(step, ir->nstorireout);
1768
1769                 print_ebin(outf->fp_ene, do_ene, do_dr, do_or, do_log ? fplog : NULL,
1770                            step, t,
1771                            eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1772             }
1773             if (ir->ePull != epullNO)
1774             {
1775                 pull_print_output(ir->pull, step, t);
1776             }
1777
1778             if (do_per_step(step, ir->nstlog))
1779             {
1780                 if (fflush(fplog) != 0)
1781                 {
1782                     gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1783                 }
1784             }
1785         }
1786         if (bDoExpanded)
1787         {
1788             /* Have to do this part _after_ outputting the logfile and the edr file */
1789             /* Gets written into the state at the beginning of next loop*/
1790             state->fep_state = lamnew;
1791         }
1792
1793         /* Remaining runtime */
1794         if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1795         {
1796             if (shellfc)
1797             {
1798                 fprintf(stderr, "\n");
1799             }
1800             print_time(stderr, runtime, step, ir, cr);
1801         }
1802
1803         /* Replica exchange */
1804         bExchanged = FALSE;
1805         if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1806             do_per_step(step, repl_ex_nst))
1807         {
1808             bExchanged = replica_exchange(fplog, cr, repl_ex,
1809                                           state_global, enerd,
1810                                           state, step, t);
1811
1812             if (bExchanged && DOMAINDECOMP(cr))
1813             {
1814                 dd_partition_system(fplog, step, cr, TRUE, 1,
1815                                     state_global, top_global, ir,
1816                                     state, &f, mdatoms, top, fr,
1817                                     vsite, shellfc, constr,
1818                                     nrnb, wcycle, FALSE);
1819             }
1820         }
1821
1822         bFirstStep       = FALSE;
1823         bInitStep        = FALSE;
1824         bStartingFromCpt = FALSE;
1825
1826         /* #######  SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1827         /* With all integrators, except VV, we need to retain the pressure
1828          * at the current step for coupling at the next step.
1829          */
1830         if ((state->flags & (1<<estPRES_PREV)) &&
1831             (bGStatEveryStep ||
1832              (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1833         {
1834             /* Store the pressure in t_state for pressure coupling
1835              * at the next MD step.
1836              */
1837             copy_mat(pres, state->pres_prev);
1838         }
1839
1840         /* #######  END SET VARIABLES FOR NEXT ITERATION ###### */
1841
1842         if ( (membed != NULL) && (!bLastStep) )
1843         {
1844             rescale_membed(step_rel, membed, state_global->x);
1845         }
1846
1847         if (bRerunMD)
1848         {
1849             if (MASTER(cr))
1850             {
1851                 /* read next frame from input trajectory */
1852                 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1853             }
1854
1855             if (PAR(cr))
1856             {
1857                 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1858             }
1859         }
1860
1861         if (!bRerunMD || !rerun_fr.bStep)
1862         {
1863             /* increase the MD step number */
1864             step++;
1865             step_rel++;
1866         }
1867
1868         cycles = wallcycle_stop(wcycle, ewcSTEP);
1869         if (DOMAINDECOMP(cr) && wcycle)
1870         {
1871             dd_cycles_add(cr->dd, cycles, ddCyclStep);
1872         }
1873
1874         if (bPMETuneRunning || bPMETuneTry)
1875         {
1876             /* PME grid + cut-off optimization with GPUs or PME nodes */
1877
1878             /* Count the total cycles over the last steps */
1879             cycles_pmes += cycles;
1880
1881             /* We can only switch cut-off at NS steps */
1882             if (step % ir->nstlist == 0)
1883             {
1884                 /* PME grid + cut-off optimization with GPUs or PME nodes */
1885                 if (bPMETuneTry)
1886                 {
1887                     if (DDMASTER(cr->dd))
1888                     {
1889                         /* PME node load is too high, start tuning */
1890                         bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
1891                     }
1892                     dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
1893
1894                     if (bPMETuneRunning || step_rel > ir->nstlist*50)
1895                     {
1896                         bPMETuneTry     = FALSE;
1897                     }
1898                 }
1899                 if (bPMETuneRunning)
1900                 {
1901                     /* init_step might not be a multiple of nstlist,
1902                      * but the first cycle is always skipped anyhow.
1903                      */
1904                     bPMETuneRunning =
1905                         pme_load_balance(pme_loadbal, cr,
1906                                          (bVerbose && MASTER(cr)) ? stderr : NULL,
1907                                          fplog,
1908                                          ir, state, cycles_pmes,
1909                                          fr->ic, fr->nbv, &fr->pmedata,
1910                                          step);
1911
1912                     /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1913                     fr->ewaldcoeff = fr->ic->ewaldcoeff;
1914                     fr->rlist      = fr->ic->rlist;
1915                     fr->rlistlong  = fr->ic->rlistlong;
1916                     fr->rcoulomb   = fr->ic->rcoulomb;
1917                     fr->rvdw       = fr->ic->rvdw;
1918                 }
1919                 cycles_pmes = 0;
1920             }
1921         }
1922
1923         if (step_rel == wcycle_get_reset_counters(wcycle) ||
1924             gs.set[eglsRESETCOUNTERS] != 0)
1925         {
1926             /* Reset all the counters related to performance over the run */
1927             reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, runtime,
1928                                fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
1929             wcycle_set_reset_counters(wcycle, -1);
1930             if (!(cr->duty & DUTY_PME))
1931             {
1932                 /* Tell our PME node to reset its counters */
1933                 gmx_pme_send_resetcounters(cr, step);
1934             }
1935             /* Correct max_hours for the elapsed time */
1936             max_hours                -= run_time/(60.0*60.0);
1937             bResetCountersHalfMaxH    = FALSE;
1938             gs.set[eglsRESETCOUNTERS] = 0;
1939         }
1940
1941     }
1942     /* End of main MD loop */
1943     debug_gmx();
1944
1945     /* Stop the time */
1946     runtime_end(runtime);
1947
1948     if (bRerunMD && MASTER(cr))
1949     {
1950         close_trj(status);
1951     }
1952
1953     if (!(cr->duty & DUTY_PME))
1954     {
1955         /* Tell the PME only node to finish */
1956         gmx_pme_send_finish(cr);
1957     }
1958
1959     if (MASTER(cr))
1960     {
1961         if (ir->nstcalcenergy > 0 && !bRerunMD)
1962         {
1963             print_ebin(outf->fp_ene, FALSE, FALSE, FALSE, fplog, step, t,
1964                        eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1965         }
1966     }
1967
1968     done_mdoutf(outf);
1969
1970     debug_gmx();
1971
1972     if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
1973     {
1974         fprintf(fplog, "Average neighborlist lifetime: %.1f steps, std.dev.: %.1f steps\n", nlh.s1/nlh.nns, sqrt(nlh.s2/nlh.nns - sqr(nlh.s1/nlh.nns)));
1975         fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
1976     }
1977
1978     if (pme_loadbal != NULL)
1979     {
1980         pme_loadbal_done(pme_loadbal, cr, fplog,
1981                          fr->nbv != NULL && fr->nbv->bUseGPU);
1982     }
1983
1984     if (shellfc && fplog)
1985     {
1986         fprintf(fplog, "Fraction of iterations that converged:           %.2f %%\n",
1987                 (nconverged*100.0)/step_rel);
1988         fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1989                 tcount/step_rel);
1990     }
1991
1992     if (repl_ex_nst > 0 && MASTER(cr))
1993     {
1994         print_replica_exchange_statistics(fplog, repl_ex);
1995     }
1996
1997     runtime->nsteps_done = step_rel;
1998
1999     return 0;
2000 }