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