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