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