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