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