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