787b704910825ce58441d9c59f076b41526fce24
[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, bSumEkinhOld, bDoReplEx, 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,
328                                  top_global, n_flexible_constraints(constr),
329                                  (ir->bContinuation ||
330                                   (DOMAINDECOMP(cr) && !MASTER(cr))) ?
331                                  NULL : state_global->x);
332     if (shellfc && ir->nstcalcenergy != 1)
333     {
334         gmx_fatal(FARGS, "You have nstcalcenergy set to a value (%d) that is different from 1.\nThis is not supported in combinations with shell particles.\nPlease make a new tpr file.", ir->nstcalcenergy);
335     }
336     if (shellfc && DOMAINDECOMP(cr))
337     {
338         gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
339     }
340     if (shellfc && ir->eI == eiNM)
341     {
342         /* Currently shells don't work with Normal Modes */
343         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");
344     }
345
346     if (vsite && ir->eI == eiNM)
347     {
348         /* Currently virtual sites don't work with Normal Modes */
349         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");
350     }
351
352     if (DEFORM(*ir))
353     {
354         tMPI_Thread_mutex_lock(&deform_init_box_mutex);
355         set_deform_reference_box(upd,
356                                  deform_init_init_step_tpx,
357                                  deform_init_box_tpx);
358         tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
359     }
360
361     {
362         double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
363         if ((io > 2000) && MASTER(cr))
364         {
365             fprintf(stderr,
366                     "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
367                     io);
368         }
369     }
370
371     if (DOMAINDECOMP(cr))
372     {
373         top = dd_init_local_top(top_global);
374
375         snew(state, 1);
376         dd_init_local_state(cr->dd, state_global, state);
377
378         if (DDMASTER(cr->dd) && ir->nstfout)
379         {
380             snew(f_global, state_global->natoms);
381         }
382     }
383     else
384     {
385         top = gmx_mtop_generate_local_top(top_global, ir);
386
387         forcerec_set_excl_load(fr, top);
388
389         state    = serial_init_local_state(state_global);
390         f_global = f;
391
392         atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
393
394         if (vsite)
395         {
396             set_vsite_top(vsite, top, mdatoms, cr);
397         }
398
399         if (ir->ePBC != epbcNONE && !fr->bMolPBC)
400         {
401             graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
402         }
403
404         if (shellfc)
405         {
406             make_local_shells(cr, mdatoms, shellfc);
407         }
408
409         setup_bonded_threading(fr, &top->idef);
410     }
411
412     /* Set up interactive MD (IMD) */
413     init_IMD(ir, cr, top_global, fplog, ir->nstcalcenergy, state_global->x,
414              nfile, fnm, oenv, imdport, Flags);
415
416     if (DOMAINDECOMP(cr))
417     {
418         /* Distribute the charge groups over the nodes from the master node */
419         dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
420                             state_global, top_global, ir,
421                             state, &f, mdatoms, top, fr,
422                             vsite, shellfc, constr,
423                             nrnb, wcycle, FALSE);
424
425     }
426
427     update_mdatoms(mdatoms, state->lambda[efptMASS]);
428
429     if (opt2bSet("-cpi", nfile, fnm))
430     {
431         bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
432     }
433     else
434     {
435         bStateFromCP = FALSE;
436     }
437
438     if (ir->bExpanded)
439     {
440         init_expanded_ensemble(bStateFromCP, ir, &state->dfhist);
441     }
442
443     if (MASTER(cr))
444     {
445         if (bStateFromCP)
446         {
447             /* Update mdebin with energy history if appending to output files */
448             if (Flags & MD_APPENDFILES)
449             {
450                 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
451             }
452             else
453             {
454                 /* We might have read an energy history from checkpoint,
455                  * free the allocated memory and reset the counts.
456                  */
457                 done_energyhistory(&state_global->enerhist);
458                 init_energyhistory(&state_global->enerhist);
459             }
460         }
461         /* Set the initial energy history in state by updating once */
462         update_energyhistory(&state_global->enerhist, mdebin);
463     }
464
465     /* Initialize constraints */
466     if (constr && !DOMAINDECOMP(cr))
467     {
468         set_constraints(constr, top, ir, mdatoms, cr);
469     }
470
471     if (repl_ex_nst > 0 && MASTER(cr))
472     {
473         repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
474                                         repl_ex_nst, repl_ex_nex, repl_ex_seed);
475     }
476
477     /* PME tuning is only supported with GPUs or PME nodes and not with rerun.
478      * PME tuning is not supported with PME only for LJ and not for Coulomb.
479      */
480     if ((Flags & MD_TUNEPME) &&
481         EEL_PME(fr->eeltype) &&
482         ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
483         !bRerunMD)
484     {
485         pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
486         cycles_pmes = 0;
487         if (cr->duty & DUTY_PME)
488         {
489             /* Start tuning right away, as we can't measure the load */
490             bPMETuneRunning = TRUE;
491         }
492         else
493         {
494             /* Separate PME nodes, we can measure the PP/PME load balance */
495             bPMETuneTry = TRUE;
496         }
497     }
498
499     if (!ir->bContinuation && !bRerunMD)
500     {
501         if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
502         {
503             /* Set the velocities of frozen particles to zero */
504             for (i = 0; i < mdatoms->homenr; i++)
505             {
506                 for (m = 0; m < DIM; m++)
507                 {
508                     if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
509                     {
510                         state->v[i][m] = 0;
511                     }
512                 }
513             }
514         }
515
516         if (constr)
517         {
518             /* Constrain the initial coordinates and velocities */
519             do_constrain_first(fplog, constr, ir, mdatoms, state,
520                                cr, nrnb, fr, top);
521         }
522         if (vsite)
523         {
524             /* Construct the virtual sites for the initial configuration */
525             construct_vsites(vsite, state->x, ir->delta_t, NULL,
526                              top->idef.iparams, top->idef.il,
527                              fr->ePBC, fr->bMolPBC, cr, state->box);
528         }
529     }
530
531     debug_gmx();
532
533     /* set free energy calculation frequency as the minimum
534        greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
535     nstfep = ir->fepvals->nstdhdl;
536     if (ir->bExpanded)
537     {
538         nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
539     }
540     if (repl_ex_nst > 0)
541     {
542         nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
543     }
544
545     /* I'm assuming we need global communication the first time! MRS */
546     cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
547                   | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
548                   | (bVV ? CGLO_PRESSURE : 0)
549                   | (bVV ? CGLO_CONSTRAINT : 0)
550                   | (bRerunMD ? CGLO_RERUNMD : 0)
551                   | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
552
553     bSumEkinhOld = FALSE;
554     compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
555                     NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
556                     constr, NULL, FALSE, state->box,
557                     top_global, &bSumEkinhOld, cglo_flags);
558     if (ir->eI == eiVVAK)
559     {
560         /* a second call to get the half step temperature initialized as well */
561         /* we do the same call as above, but turn the pressure off -- internally to
562            compute_globals, this is recognized as a velocity verlet half-step
563            kinetic energy calculation.  This minimized excess variables, but
564            perhaps loses some logic?*/
565
566         compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
567                         NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
568                         constr, NULL, FALSE, state->box,
569                         top_global, &bSumEkinhOld,
570                         cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
571     }
572
573     /* Calculate the initial half step temperature, and save the ekinh_old */
574     if (!(Flags & MD_STARTFROMCPT))
575     {
576         for (i = 0; (i < ir->opts.ngtc); i++)
577         {
578             copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
579         }
580     }
581     if (ir->eI != eiVV)
582     {
583         enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
584                                      and there is no previous step */
585     }
586
587     /* if using an iterative algorithm, we need to create a working directory for the state. */
588     if (bIterativeCase)
589     {
590         bufstate = init_bufstate(state);
591     }
592
593     /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
594        temperature control */
595     trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
596
597     if (MASTER(cr))
598     {
599         if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
600         {
601             fprintf(fplog,
602                     "RMS relative constraint deviation after constraining: %.2e\n",
603                     constr_rmsd(constr, FALSE));
604         }
605         if (EI_STATE_VELOCITY(ir->eI))
606         {
607             fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
608         }
609         if (bRerunMD)
610         {
611             fprintf(stderr, "starting md rerun '%s', reading coordinates from"
612                     " input trajectory '%s'\n\n",
613                     *(top_global->name), opt2fn("-rerun", nfile, fnm));
614             if (bVerbose)
615             {
616                 fprintf(stderr, "Calculated time to finish depends on nsteps from "
617                         "run input file,\nwhich may not correspond to the time "
618                         "needed to process input trajectory.\n\n");
619             }
620         }
621         else
622         {
623             char tbuf[20];
624             fprintf(stderr, "starting mdrun '%s'\n",
625                     *(top_global->name));
626             if (ir->nsteps >= 0)
627             {
628                 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
629             }
630             else
631             {
632                 sprintf(tbuf, "%s", "infinite");
633             }
634             if (ir->init_step > 0)
635             {
636                 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
637                         gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
638                         gmx_step_str(ir->init_step, sbuf2),
639                         ir->init_step*ir->delta_t);
640             }
641             else
642             {
643                 fprintf(stderr, "%s steps, %s ps.\n",
644                         gmx_step_str(ir->nsteps, sbuf), tbuf);
645             }
646         }
647         fprintf(fplog, "\n");
648     }
649
650     walltime_accounting_start(walltime_accounting);
651     wallcycle_start(wcycle, ewcRUN);
652     print_start(fplog, cr, walltime_accounting, "mdrun");
653
654     /* safest point to do file checkpointing is here.  More general point would be immediately before integrator call */
655 #ifdef GMX_FAHCORE
656     chkpt_ret = fcCheckPointParallel( cr->nodeid,
657                                       NULL, 0);
658     if (chkpt_ret == 0)
659     {
660         gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
661     }
662 #endif
663
664     debug_gmx();
665     /***********************************************************
666      *
667      *             Loop over MD steps
668      *
669      ************************************************************/
670
671     /* if rerunMD then read coordinates and velocities from input trajectory */
672     if (bRerunMD)
673     {
674         if (getenv("GMX_FORCE_UPDATE"))
675         {
676             bForceUpdate = TRUE;
677         }
678
679         rerun_fr.natoms = 0;
680         if (MASTER(cr))
681         {
682             bNotLastFrame = read_first_frame(oenv, &status,
683                                              opt2fn("-rerun", nfile, fnm),
684                                              &rerun_fr, TRX_NEED_X | TRX_READ_V);
685             if (rerun_fr.natoms != top_global->natoms)
686             {
687                 gmx_fatal(FARGS,
688                           "Number of atoms in trajectory (%d) does not match the "
689                           "run input file (%d)\n",
690                           rerun_fr.natoms, top_global->natoms);
691             }
692             if (ir->ePBC != epbcNONE)
693             {
694                 if (!rerun_fr.bBox)
695                 {
696                     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);
697                 }
698                 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
699                 {
700                     gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
701                 }
702             }
703         }
704
705         if (PAR(cr))
706         {
707             rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
708         }
709
710         if (ir->ePBC != epbcNONE)
711         {
712             /* Set the shift vectors.
713              * Necessary here when have a static box different from the tpr box.
714              */
715             calc_shifts(rerun_fr.box, fr->shift_vec);
716         }
717     }
718
719     /* loop over MD steps or if rerunMD to end of input trajectory */
720     bFirstStep = TRUE;
721     /* Skip the first Nose-Hoover integration when we get the state from tpx */
722     bStateFromTPX    = !bStateFromCP;
723     bInitStep        = bFirstStep && (bStateFromTPX || bVV);
724     bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
725     bSumEkinhOld     = FALSE;
726     bExchanged       = FALSE;
727     bNeedRepartition = FALSE;
728
729     init_global_signals(&gs, cr, ir, repl_ex_nst);
730
731     step     = ir->init_step;
732     step_rel = 0;
733
734     init_nlistheuristics(&nlh, bGStatEveryStep, step);
735
736     if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
737     {
738         /* check how many steps are left in other sims */
739         multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
740     }
741
742
743     /* and stop now if we should */
744     bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
745                  ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
746     while (!bLastStep || (bRerunMD && bNotLastFrame))
747     {
748
749         wallcycle_start(wcycle, ewcSTEP);
750
751         if (bRerunMD)
752         {
753             if (rerun_fr.bStep)
754             {
755                 step     = rerun_fr.step;
756                 step_rel = step - ir->init_step;
757             }
758             if (rerun_fr.bTime)
759             {
760                 t = rerun_fr.time;
761             }
762             else
763             {
764                 t = step;
765             }
766         }
767         else
768         {
769             bLastStep = (step_rel == ir->nsteps);
770             t         = t0 + step*ir->delta_t;
771         }
772
773         if (ir->efep != efepNO || ir->bSimTemp)
774         {
775             /* find and set the current lambdas.  If rerunning, we either read in a state, or a lambda value,
776                requiring different logic. */
777
778             set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
779             bDoDHDL      = do_per_step(step, ir->fepvals->nstdhdl);
780             bDoFEP       = (do_per_step(step, nstfep) && (ir->efep != efepNO));
781             bDoExpanded  = (do_per_step(step, ir->expandedvals->nstexpanded)
782                             && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
783         }
784
785         bDoReplEx = ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
786                      do_per_step(step, repl_ex_nst));
787
788         if (bSimAnn)
789         {
790             update_annealing_target_temp(&(ir->opts), t);
791         }
792
793         if (bRerunMD)
794         {
795             if (!DOMAINDECOMP(cr) || MASTER(cr))
796             {
797                 for (i = 0; i < state_global->natoms; i++)
798                 {
799                     copy_rvec(rerun_fr.x[i], state_global->x[i]);
800                 }
801                 if (rerun_fr.bV)
802                 {
803                     for (i = 0; i < state_global->natoms; i++)
804                     {
805                         copy_rvec(rerun_fr.v[i], state_global->v[i]);
806                     }
807                 }
808                 else
809                 {
810                     for (i = 0; i < state_global->natoms; i++)
811                     {
812                         clear_rvec(state_global->v[i]);
813                     }
814                     if (bRerunWarnNoV)
815                     {
816                         fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
817                                 "         Ekin, temperature and pressure are incorrect,\n"
818                                 "         the virial will be incorrect when constraints are present.\n"
819                                 "\n");
820                         bRerunWarnNoV = FALSE;
821                     }
822                 }
823             }
824             copy_mat(rerun_fr.box, state_global->box);
825             copy_mat(state_global->box, state->box);
826
827             if (vsite && (Flags & MD_RERUN_VSITE))
828             {
829                 if (DOMAINDECOMP(cr))
830                 {
831                     gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
832                 }
833                 if (graph)
834                 {
835                     /* Following is necessary because the graph may get out of sync
836                      * with the coordinates if we only have every N'th coordinate set
837                      */
838                     mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
839                     shift_self(graph, state->box, state->x);
840                 }
841                 construct_vsites(vsite, state->x, ir->delta_t, state->v,
842                                  top->idef.iparams, top->idef.il,
843                                  fr->ePBC, fr->bMolPBC, cr, state->box);
844                 if (graph)
845                 {
846                     unshift_self(graph, state->box, state->x);
847                 }
848             }
849         }
850
851         /* Stop Center of Mass motion */
852         bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
853
854         if (bRerunMD)
855         {
856             /* for rerun MD always do Neighbour Searching */
857             bNS      = (bFirstStep || ir->nstlist != 0);
858             bNStList = bNS;
859         }
860         else
861         {
862             /* Determine whether or not to do Neighbour Searching and LR */
863             bNStList = (ir->nstlist > 0  && step % ir->nstlist == 0);
864
865             bNS = (bFirstStep || bExchanged || bNeedRepartition || bNStList || bDoFEP ||
866                    (ir->nstlist == -1 && nlh.nabnsb > 0));
867
868             if (bNS && ir->nstlist == -1)
869             {
870                 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bNeedRepartition || bDoFEP, step);
871             }
872         }
873
874         /* check whether we should stop because another simulation has
875            stopped. */
876         if (MULTISIM(cr))
877         {
878             if ( (multisim_nsteps >= 0) &&  (step_rel >= multisim_nsteps)  &&
879                  (multisim_nsteps != ir->nsteps) )
880             {
881                 if (bNS)
882                 {
883                     if (MASTER(cr))
884                     {
885                         fprintf(stderr,
886                                 "Stopping simulation %d because another one has finished\n",
887                                 cr->ms->sim);
888                     }
889                     bLastStep         = TRUE;
890                     gs.sig[eglsCHKPT] = 1;
891                 }
892             }
893         }
894
895         /* < 0 means stop at next step, > 0 means stop at next NS step */
896         if ( (gs.set[eglsSTOPCOND] < 0) ||
897              ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
898         {
899             bLastStep = TRUE;
900         }
901
902         /* Determine whether or not to update the Born radii if doing GB */
903         bBornRadii = bFirstStep;
904         if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
905         {
906             bBornRadii = TRUE;
907         }
908
909         do_log     = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
910         do_verbose = bVerbose &&
911             (step % stepout == 0 || bFirstStep || bLastStep);
912
913         if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
914         {
915             if (bRerunMD)
916             {
917                 bMasterState = TRUE;
918             }
919             else
920             {
921                 bMasterState = FALSE;
922                 /* Correct the new box if it is too skewed */
923                 if (DYNAMIC_BOX(*ir))
924                 {
925                     if (correct_box(fplog, step, state->box, graph))
926                     {
927                         bMasterState = TRUE;
928                     }
929                 }
930                 if (DOMAINDECOMP(cr) && bMasterState)
931                 {
932                     dd_collect_state(cr->dd, state, state_global);
933                 }
934             }
935
936             if (DOMAINDECOMP(cr))
937             {
938                 /* Repartition the domain decomposition */
939                 wallcycle_start(wcycle, ewcDOMDEC);
940                 dd_partition_system(fplog, step, cr,
941                                     bMasterState, nstglobalcomm,
942                                     state_global, top_global, ir,
943                                     state, &f, mdatoms, top, fr,
944                                     vsite, shellfc, constr,
945                                     nrnb, wcycle,
946                                     do_verbose && !bPMETuneRunning);
947                 wallcycle_stop(wcycle, ewcDOMDEC);
948                 /* If using an iterative integrator, reallocate space to match the decomposition */
949             }
950         }
951
952         if (MASTER(cr) && do_log)
953         {
954             print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
955         }
956
957         if (ir->efep != efepNO)
958         {
959             update_mdatoms(mdatoms, state->lambda[efptMASS]);
960         }
961
962         if ((bRerunMD && rerun_fr.bV) || bExchanged)
963         {
964
965             /* We need the kinetic energy at minus the half step for determining
966              * the full step kinetic energy and possibly for T-coupling.*/
967             /* This may not be quite working correctly yet . . . . */
968             compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
969                             wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
970                             constr, NULL, FALSE, state->box,
971                             top_global, &bSumEkinhOld,
972                             CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
973         }
974         clear_mat(force_vir);
975
976         /* We write a checkpoint at this MD step when:
977          * either at an NS step when we signalled through gs,
978          * or at the last step (but not when we do not want confout),
979          * but never at the first step or with rerun.
980          */
981         bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
982                  (bLastStep && (Flags & MD_CONFOUT))) &&
983                 step > ir->init_step && !bRerunMD);
984         if (bCPT)
985         {
986             gs.set[eglsCHKPT] = 0;
987         }
988
989         /* Determine the energy and pressure:
990          * at nstcalcenergy steps and at energy output steps (set below).
991          */
992         if (EI_VV(ir->eI) && (!bInitStep))
993         {
994             /* for vv, the first half of the integration actually corresponds
995                to the previous step.  bCalcEner is only required to be evaluated on the 'next' step,
996                but the virial needs to be calculated on both the current step and the 'next' step. Future
997                reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
998
999             bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1000             bCalcVir  = bCalcEner ||
1001                 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1002         }
1003         else
1004         {
1005             bCalcEner = do_per_step(step, ir->nstcalcenergy);
1006             bCalcVir  = bCalcEner ||
1007                 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1008         }
1009
1010         /* Do we need global communication ? */
1011         bGStat = (bCalcVir || bCalcEner || bStopCM ||
1012                   do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1013                   (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1014
1015         do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1016
1017         if (do_ene || do_log || bDoReplEx)
1018         {
1019             bCalcVir  = TRUE;
1020             bCalcEner = TRUE;
1021             bGStat    = TRUE;
1022         }
1023
1024         /* these CGLO_ options remain the same throughout the iteration */
1025         cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1026                       (bGStat ? CGLO_GSTAT : 0)
1027                       );
1028
1029         force_flags = (GMX_FORCE_STATECHANGED |
1030                        ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1031                        GMX_FORCE_ALLFORCES |
1032                        GMX_FORCE_SEPLRF |
1033                        (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1034                        (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1035                        (bDoFEP ? GMX_FORCE_DHDL : 0)
1036                        );
1037
1038         if (fr->bTwinRange)
1039         {
1040             if (do_per_step(step, ir->nstcalclr))
1041             {
1042                 force_flags |= GMX_FORCE_DO_LR;
1043             }
1044         }
1045
1046         if (shellfc)
1047         {
1048             /* Now is the time to relax the shells */
1049             count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1050                                         ir, bNS, force_flags,
1051                                         top,
1052                                         constr, enerd, fcd,
1053                                         state, f, force_vir, mdatoms,
1054                                         nrnb, wcycle, graph, groups,
1055                                         shellfc, fr, bBornRadii, t, mu_tot,
1056                                         &bConverged, vsite,
1057                                         mdoutf_get_fp_field(outf));
1058             tcount += count;
1059
1060             if (bConverged)
1061             {
1062                 nconverged++;
1063             }
1064         }
1065         else
1066         {
1067             /* The coordinates (x) are shifted (to get whole molecules)
1068              * in do_force.
1069              * This is parallellized as well, and does communication too.
1070              * Check comments in sim_util.c
1071              */
1072             do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1073                      state->box, state->x, &state->hist,
1074                      f, force_vir, mdatoms, enerd, fcd,
1075                      state->lambda, graph,
1076                      fr, vsite, mu_tot, t, mdoutf_get_fp_field(outf), ed, bBornRadii,
1077                      (bNS ? GMX_FORCE_NS : 0) | force_flags);
1078         }
1079
1080         if (bVV && !bStartingFromCpt && !bRerunMD)
1081         /*  ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1082         {
1083             if (ir->eI == eiVV && bInitStep)
1084             {
1085                 /* if using velocity verlet with full time step Ekin,
1086                  * take the first half step only to compute the
1087                  * virial for the first step. From there,
1088                  * revert back to the initial coordinates
1089                  * so that the input is actually the initial step.
1090                  */
1091                 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1092             }
1093             else
1094             {
1095                 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1096                 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1097             }
1098
1099             /* If we are using twin-range interactions where the long-range component
1100              * is only evaluated every nstcalclr>1 steps, we should do a special update
1101              * step to combine the long-range forces on these steps.
1102              * For nstcalclr=1 this is not done, since the forces would have been added
1103              * directly to the short-range forces already.
1104              *
1105              * TODO Remove various aspects of VV+twin-range in master
1106              * branch, because VV integrators did not ever support
1107              * twin-range multiple time stepping with constraints.
1108              */
1109             bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1110
1111             update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1112                           f, bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1113                           ekind, M, upd, bInitStep, etrtVELOCITY1,
1114                           cr, nrnb, constr, &top->idef);
1115
1116             if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1117             {
1118                 gmx_iterate_init(&iterate, TRUE);
1119             }
1120             /* for iterations, we save these vectors, as we will be self-consistently iterating
1121                the calculations */
1122
1123             /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1124
1125             /* save the state */
1126             if (iterate.bIterationActive)
1127             {
1128                 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1129             }
1130
1131             bFirstIterate = TRUE;
1132             while (bFirstIterate || iterate.bIterationActive)
1133             {
1134                 if (iterate.bIterationActive)
1135                 {
1136                     copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1137                     if (bFirstIterate && bTrotter)
1138                     {
1139                         /* The first time through, we need a decent first estimate
1140                            of veta(t+dt) to compute the constraints.  Do
1141                            this by computing the box volume part of the
1142                            trotter integration at this time. Nothing else
1143                            should be changed by this routine here.  If
1144                            !(first time), we start with the previous value
1145                            of veta.  */
1146
1147                         veta_save = state->veta;
1148                         trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1149                         vetanew     = state->veta;
1150                         state->veta = veta_save;
1151                     }
1152                 }
1153
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 (bCalcVir && bUpdateDoLR && ir->nstcalclr > 1)
1163                     {
1164                         /* Correct the virial for multiple time stepping */
1165                         m_sub(shake_vir, fr->vir_twin_constr, shake_vir);
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             dvdl_constr = 0;
1464
1465             if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1466             {
1467                 wallcycle_start(wcycle, ewcUPDATE);
1468                 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1469                 if (bTrotter)
1470                 {
1471                     if (iterate.bIterationActive)
1472                     {
1473                         if (bFirstIterate)
1474                         {
1475                             scalevir = 1;
1476                         }
1477                         else
1478                         {
1479                             /* we use a new value of scalevir to converge the iterations faster */
1480                             scalevir = tracevir/trace(shake_vir);
1481                         }
1482                         msmul(shake_vir, scalevir, shake_vir);
1483                         m_add(force_vir, shake_vir, total_vir);
1484                         clear_mat(shake_vir);
1485                     }
1486                     trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1487                     /* We can only do Berendsen coupling after we have summed
1488                      * the kinetic energy or virial. Since the happens
1489                      * in global_state after update, we should only do it at
1490                      * step % nstlist = 1 with bGStatEveryStep=FALSE.
1491                      */
1492                 }
1493                 else
1494                 {
1495                     update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1496                     update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1497                 }
1498
1499                 if (bVV)
1500                 {
1501                     bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1502
1503                     /* velocity half-step update */
1504                     update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1505                                   bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1506                                   ekind, M, upd, FALSE, etrtVELOCITY2,
1507                                   cr, nrnb, constr, &top->idef);
1508                 }
1509
1510                 /* Above, initialize just copies ekinh into ekin,
1511                  * it doesn't copy position (for VV),
1512                  * and entire integrator for MD.
1513                  */
1514
1515                 if (ir->eI == eiVVAK)
1516                 {
1517                     copy_rvecn(state->x, cbuf, 0, state->natoms);
1518                 }
1519                 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1520
1521                 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1522                               bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1523                               ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1524                 wallcycle_stop(wcycle, ewcUPDATE);
1525
1526                 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1527                                    fr->bMolPBC, graph, f,
1528                                    &top->idef, shake_vir,
1529                                    cr, nrnb, wcycle, upd, constr,
1530                                    FALSE, bCalcVir, state->veta);
1531
1532                 if (bCalcVir && bUpdateDoLR && ir->nstcalclr > 1)
1533                 {
1534                     /* Correct the virial for multiple time stepping */
1535                     m_sub(shake_vir, fr->vir_twin_constr, shake_vir);
1536                 }
1537
1538                 if (ir->eI == eiVVAK)
1539                 {
1540                     /* erase F_EKIN and F_TEMP here? */
1541                     /* just compute the kinetic energy at the half step to perform a trotter step */
1542                     compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1543                                     wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1544                                     constr, NULL, FALSE, lastbox,
1545                                     top_global, &bSumEkinhOld,
1546                                     cglo_flags | CGLO_TEMPERATURE
1547                                     );
1548                     wallcycle_start(wcycle, ewcUPDATE);
1549                     trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1550                     /* now we know the scaling, we can compute the positions again again */
1551                     copy_rvecn(cbuf, state->x, 0, state->natoms);
1552
1553                     bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1554
1555                     update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1556                                   bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1557                                   ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1558                     wallcycle_stop(wcycle, ewcUPDATE);
1559
1560                     /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1561                     /* are the small terms in the shake_vir here due
1562                      * to numerical errors, or are they important
1563                      * physically? I'm thinking they are just errors, but not completely sure.
1564                      * For now, will call without actually constraining, constr=NULL*/
1565                     update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1566                                        state, fr->bMolPBC, graph, f,
1567                                        &top->idef, tmp_vir,
1568                                        cr, nrnb, wcycle, upd, NULL,
1569                                        FALSE, bCalcVir,
1570                                        state->veta);
1571                 }
1572
1573                 if (fr->bSepDVDL && fplog && do_log)
1574                 {
1575                     gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1576                 }
1577                 if (bVV)
1578                 {
1579                     /* this factor or 2 correction is necessary
1580                        because half of the constraint force is removed
1581                        in the vv step, so we have to double it.  See
1582                        the Redmine issue #1255.  It is not yet clear
1583                        if the factor of 2 is exact, or just a very
1584                        good approximation, and this will be
1585                        investigated.  The next step is to see if this
1586                        can be done adding a dhdl contribution from the
1587                        rattle step, but this is somewhat more
1588                        complicated with the current code. Will be
1589                        investigated, hopefully for 4.6.3. However,
1590                        this current solution is much better than
1591                        having it completely wrong.
1592                      */
1593                     enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1594                 }
1595                 else
1596                 {
1597                     enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1598                 }
1599             }
1600             else if (graph)
1601             {
1602                 /* Need to unshift here */
1603                 unshift_self(graph, state->box, state->x);
1604             }
1605
1606             if (vsite != NULL)
1607             {
1608                 wallcycle_start(wcycle, ewcVSITECONSTR);
1609                 if (graph != NULL)
1610                 {
1611                     shift_self(graph, state->box, state->x);
1612                 }
1613                 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1614                                  top->idef.iparams, top->idef.il,
1615                                  fr->ePBC, fr->bMolPBC, cr, state->box);
1616
1617                 if (graph != NULL)
1618                 {
1619                     unshift_self(graph, state->box, state->x);
1620                 }
1621                 wallcycle_stop(wcycle, ewcVSITECONSTR);
1622             }
1623
1624             /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints  ############ */
1625             /* With Leap-Frog we can skip compute_globals at
1626              * non-communication steps, but we need to calculate
1627              * the kinetic energy one step before communication.
1628              */
1629             if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1630             {
1631                 if (ir->nstlist == -1 && bFirstIterate)
1632                 {
1633                     gs.sig[eglsNABNSB] = nlh.nabnsb;
1634                 }
1635                 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1636                                 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1637                                 constr,
1638                                 bFirstIterate ? &gs : NULL,
1639                                 (step_rel % gs.nstms == 0) &&
1640                                 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1641                                 lastbox,
1642                                 top_global, &bSumEkinhOld,
1643                                 cglo_flags
1644                                 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1645                                 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1646                                 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1647                                 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1648                                 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1649                                 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1650                                 | CGLO_CONSTRAINT
1651                                 );
1652                 if (ir->nstlist == -1 && bFirstIterate)
1653                 {
1654                     nlh.nabnsb         = gs.set[eglsNABNSB];
1655                     gs.set[eglsNABNSB] = 0;
1656                 }
1657             }
1658             /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1659             /* #############  END CALC EKIN AND PRESSURE ################# */
1660
1661             /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1662                the virial that should probably be addressed eventually. state->veta has better properies,
1663                but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1664                generate the new shake_vir, but test the veta value for convergence.  This will take some thought. */
1665
1666             if (iterate.bIterationActive &&
1667                 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1668                                trace(shake_vir), &tracevir))
1669             {
1670                 break;
1671             }
1672             bFirstIterate = FALSE;
1673         }
1674
1675         if (!bVV || bRerunMD)
1676         {
1677             /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1678             sum_dhdl(enerd, state->lambda, ir->fepvals);
1679         }
1680         update_box(fplog, step, ir, mdatoms, state, f,
1681                    ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1682
1683         /* ################# END UPDATE STEP 2 ################# */
1684         /* #### We now have r(t+dt) and v(t+dt/2)  ############# */
1685
1686         /* The coordinates (x) were unshifted in update */
1687         if (!bGStat)
1688         {
1689             /* We will not sum ekinh_old,
1690              * so signal that we still have to do it.
1691              */
1692             bSumEkinhOld = TRUE;
1693         }
1694
1695         /* #########  BEGIN PREPARING EDR OUTPUT  ###########  */
1696
1697         /* use the directly determined last velocity, not actually the averaged half steps */
1698         if (bTrotter && ir->eI == eiVV)
1699         {
1700             enerd->term[F_EKIN] = last_ekin;
1701         }
1702         enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1703
1704         if (bVV)
1705         {
1706             enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1707         }
1708         else
1709         {
1710             enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1711         }
1712         /* #########  END PREPARING EDR OUTPUT  ###########  */
1713
1714         /* Output stuff */
1715         if (MASTER(cr))
1716         {
1717             gmx_bool do_dr, do_or;
1718
1719             if (fplog && do_log && bDoExpanded)
1720             {
1721                 /* only needed if doing expanded ensemble */
1722                 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1723                                           &state_global->dfhist, state->fep_state, ir->nstlog, step);
1724             }
1725             if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1726             {
1727                 if (bCalcEner)
1728                 {
1729                     upd_mdebin(mdebin, bDoDHDL, TRUE,
1730                                t, mdatoms->tmass, enerd, state,
1731                                ir->fepvals, ir->expandedvals, lastbox,
1732                                shake_vir, force_vir, total_vir, pres,
1733                                ekind, mu_tot, constr);
1734                 }
1735                 else
1736                 {
1737                     upd_mdebin_step(mdebin);
1738                 }
1739
1740                 do_dr  = do_per_step(step, ir->nstdisreout);
1741                 do_or  = do_per_step(step, ir->nstorireout);
1742
1743                 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1744                            step, t,
1745                            eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1746             }
1747             if (ir->ePull != epullNO)
1748             {
1749                 pull_print_output(ir->pull, step, t);
1750             }
1751
1752             if (do_per_step(step, ir->nstlog))
1753             {
1754                 if (fflush(fplog) != 0)
1755                 {
1756                     gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1757                 }
1758             }
1759         }
1760         if (bDoExpanded)
1761         {
1762             /* Have to do this part _after_ outputting the logfile and the edr file */
1763             /* Gets written into the state at the beginning of next loop*/
1764             state->fep_state = lamnew;
1765         }
1766         /* Print the remaining wall clock time for the run */
1767         if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1768         {
1769             if (shellfc)
1770             {
1771                 fprintf(stderr, "\n");
1772             }
1773             print_time(stderr, walltime_accounting, step, ir, cr);
1774         }
1775
1776         /* Ion/water position swapping.
1777          * Not done in last step since trajectory writing happens before this call
1778          * in the MD loop and exchanges would be lost anyway. */
1779         bNeedRepartition = FALSE;
1780         if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1781             do_per_step(step, ir->swap->nstswap))
1782         {
1783             bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1784                                              bRerunMD ? rerun_fr.x   : state->x,
1785                                              bRerunMD ? rerun_fr.box : state->box,
1786                                              top_global, MASTER(cr) && bVerbose, bRerunMD);
1787
1788             if (bNeedRepartition && DOMAINDECOMP(cr))
1789             {
1790                 dd_collect_state(cr->dd, state, state_global);
1791             }
1792         }
1793
1794         /* Replica exchange */
1795         bExchanged = FALSE;
1796         if (bDoReplEx)
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 }