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