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