Add fatal errors for VV and twin-range MTS
[alexxy/gromacs.git] / src / mdlib / md_support.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  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012,2013, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include "typedefs.h"
43 #include "string2.h"
44 #include "smalloc.h"
45 #include "mdrun.h"
46 #include "domdec.h"
47 #include "mtop_util.h"
48 #include "gmx_wallcycle.h"
49 #include "vcm.h"
50 #include "nrnb.h"
51 #include "md_logging.h"
52 #include "md_support.h"
53 #include "names.h"
54
55 /* Is the signal in one simulation independent of other simulations? */
56 gmx_bool gs_simlocal[eglsNR] = { TRUE, FALSE, FALSE, TRUE };
57
58 /* check which of the multisim simulations has the shortest number of
59    steps and return that number of nsteps */
60 gmx_large_int_t get_multisim_nsteps(const t_commrec *cr,
61                                     gmx_large_int_t  nsteps)
62 {
63     gmx_large_int_t steps_out;
64
65     if MASTER(cr)
66     {
67         gmx_large_int_t *buf;
68         int              s;
69
70         snew(buf, cr->ms->nsim);
71
72         buf[cr->ms->sim] = nsteps;
73         gmx_sumli_sim(cr->ms->nsim, buf, cr->ms);
74
75         steps_out = -1;
76         for (s = 0; s < cr->ms->nsim; s++)
77         {
78             /* find the smallest positive number */
79             if (buf[s] >= 0 && ((steps_out < 0) || (buf[s] < steps_out)) )
80             {
81                 steps_out = buf[s];
82             }
83         }
84         sfree(buf);
85
86         /* if we're the limiting simulation, don't do anything */
87         if (steps_out >= 0 && steps_out < nsteps)
88         {
89             char strbuf[255];
90             snprintf(strbuf, 255, "Will stop simulation %%d after %s steps (another simulation will end then).\n", gmx_large_int_pfmt);
91             fprintf(stderr, strbuf, cr->ms->sim, steps_out);
92         }
93     }
94     /* broadcast to non-masters */
95     gmx_bcast(sizeof(gmx_large_int_t), &steps_out, cr);
96     return steps_out;
97 }
98
99 int multisim_min(const gmx_multisim_t *ms, int nmin, int n)
100 {
101     int     *buf;
102     gmx_bool bPos, bEqual;
103     int      s, d;
104
105     snew(buf, ms->nsim);
106     buf[ms->sim] = n;
107     gmx_sumi_sim(ms->nsim, buf, ms);
108     bPos   = TRUE;
109     bEqual = TRUE;
110     for (s = 0; s < ms->nsim; s++)
111     {
112         bPos   = bPos   && (buf[s] > 0);
113         bEqual = bEqual && (buf[s] == buf[0]);
114     }
115     if (bPos)
116     {
117         if (bEqual)
118         {
119             nmin = min(nmin, buf[0]);
120         }
121         else
122         {
123             /* Find the least common multiple */
124             for (d = 2; d < nmin; d++)
125             {
126                 s = 0;
127                 while (s < ms->nsim && d % buf[s] == 0)
128                 {
129                     s++;
130                 }
131                 if (s == ms->nsim)
132                 {
133                     /* We found the LCM and it is less than nmin */
134                     nmin = d;
135                     break;
136                 }
137             }
138         }
139     }
140     sfree(buf);
141
142     return nmin;
143 }
144
145 int multisim_nstsimsync(const t_commrec *cr,
146                         const t_inputrec *ir, int repl_ex_nst)
147 {
148     int nmin;
149
150     if (MASTER(cr))
151     {
152         nmin = INT_MAX;
153         nmin = multisim_min(cr->ms, nmin, ir->nstlist);
154         nmin = multisim_min(cr->ms, nmin, ir->nstcalcenergy);
155         nmin = multisim_min(cr->ms, nmin, repl_ex_nst);
156         if (nmin == INT_MAX)
157         {
158             gmx_fatal(FARGS, "Can not find an appropriate interval for inter-simulation communication, since nstlist, nstcalcenergy and -replex are all <= 0");
159         }
160         /* Avoid inter-simulation communication at every (second) step */
161         if (nmin <= 2)
162         {
163             nmin = 10;
164         }
165     }
166
167     gmx_bcast(sizeof(int), &nmin, cr);
168
169     return nmin;
170 }
171
172 void init_global_signals(globsig_t *gs, const t_commrec *cr,
173                          const t_inputrec *ir, int repl_ex_nst)
174 {
175     int i;
176
177     if (MULTISIM(cr))
178     {
179         gs->nstms = multisim_nstsimsync(cr, ir, repl_ex_nst);
180         if (debug)
181         {
182             fprintf(debug, "Syncing simulations for checkpointing and termination every %d steps\n", gs->nstms);
183         }
184     }
185     else
186     {
187         gs->nstms = 1;
188     }
189
190     for (i = 0; i < eglsNR; i++)
191     {
192         gs->sig[i] = 0;
193         gs->set[i] = 0;
194     }
195 }
196
197 void copy_coupling_state(t_state *statea, t_state *stateb,
198                          gmx_ekindata_t *ekinda, gmx_ekindata_t *ekindb, t_grpopts* opts)
199 {
200
201     /* MRS note -- might be able to get rid of some of the arguments.  Look over it when it's all debugged */
202
203     int i, j, nc;
204
205     /* Make sure we have enough space for x and v */
206     if (statea->nalloc > stateb->nalloc)
207     {
208         stateb->nalloc = statea->nalloc;
209         srenew(stateb->x, stateb->nalloc);
210         srenew(stateb->v, stateb->nalloc);
211     }
212
213     stateb->natoms     = statea->natoms;
214     stateb->ngtc       = statea->ngtc;
215     stateb->nnhpres    = statea->nnhpres;
216     stateb->veta       = statea->veta;
217     if (ekinda)
218     {
219         copy_mat(ekinda->ekin, ekindb->ekin);
220         for (i = 0; i < stateb->ngtc; i++)
221         {
222             ekindb->tcstat[i].T  = ekinda->tcstat[i].T;
223             ekindb->tcstat[i].Th = ekinda->tcstat[i].Th;
224             copy_mat(ekinda->tcstat[i].ekinh, ekindb->tcstat[i].ekinh);
225             copy_mat(ekinda->tcstat[i].ekinf, ekindb->tcstat[i].ekinf);
226             ekindb->tcstat[i].ekinscalef_nhc =  ekinda->tcstat[i].ekinscalef_nhc;
227             ekindb->tcstat[i].ekinscaleh_nhc =  ekinda->tcstat[i].ekinscaleh_nhc;
228             ekindb->tcstat[i].vscale_nhc     =  ekinda->tcstat[i].vscale_nhc;
229         }
230     }
231     copy_rvecn(statea->x, stateb->x, 0, stateb->natoms);
232     copy_rvecn(statea->v, stateb->v, 0, stateb->natoms);
233     copy_mat(statea->box, stateb->box);
234     copy_mat(statea->box_rel, stateb->box_rel);
235     copy_mat(statea->boxv, stateb->boxv);
236
237     for (i = 0; i < stateb->ngtc; i++)
238     {
239         nc = i*opts->nhchainlength;
240         for (j = 0; j < opts->nhchainlength; j++)
241         {
242             stateb->nosehoover_xi[nc+j]  = statea->nosehoover_xi[nc+j];
243             stateb->nosehoover_vxi[nc+j] = statea->nosehoover_vxi[nc+j];
244         }
245     }
246     if (stateb->nhpres_xi != NULL)
247     {
248         for (i = 0; i < stateb->nnhpres; i++)
249         {
250             nc = i*opts->nhchainlength;
251             for (j = 0; j < opts->nhchainlength; j++)
252             {
253                 stateb->nhpres_xi[nc+j]  = statea->nhpres_xi[nc+j];
254                 stateb->nhpres_vxi[nc+j] = statea->nhpres_vxi[nc+j];
255             }
256         }
257     }
258 }
259
260 real compute_conserved_from_auxiliary(t_inputrec *ir, t_state *state, t_extmass *MassQ)
261 {
262     real quantity = 0;
263     switch (ir->etc)
264     {
265         case etcNO:
266             break;
267         case etcBERENDSEN:
268             break;
269         case etcNOSEHOOVER:
270             quantity = NPT_energy(ir, state, MassQ);
271             break;
272         case etcVRESCALE:
273             quantity = vrescale_energy(&(ir->opts), state->therm_integral);
274             break;
275         default:
276             break;
277     }
278     return quantity;
279 }
280
281 void compute_globals(FILE *fplog, gmx_global_stat_t gstat, t_commrec *cr, t_inputrec *ir,
282                      t_forcerec *fr, gmx_ekindata_t *ekind,
283                      t_state *state, t_state *state_global, t_mdatoms *mdatoms,
284                      t_nrnb *nrnb, t_vcm *vcm, gmx_wallcycle_t wcycle,
285                      gmx_enerdata_t *enerd, tensor force_vir, tensor shake_vir, tensor total_vir,
286                      tensor pres, rvec mu_tot, gmx_constr_t constr,
287                      globsig_t *gs, gmx_bool bInterSimGS,
288                      matrix box, gmx_mtop_t *top_global, real *pcurr,
289                      int natoms, gmx_bool *bSumEkinhOld, int flags)
290 {
291     int      i, gsi;
292     real     gs_buf[eglsNR];
293     tensor   corr_vir, corr_pres, shakeall_vir;
294     gmx_bool bEner, bPres, bTemp, bVV;
295     gmx_bool bRerunMD, bStopCM, bGStat, bIterate,
296              bFirstIterate, bReadEkin, bEkinAveVel, bScaleEkin, bConstrain;
297     real     ekin, temp, prescorr, enercorr, dvdlcorr, dvdl_ekin;
298
299     /* translate CGLO flags to gmx_booleans */
300     bRerunMD = flags & CGLO_RERUNMD;
301     bStopCM  = flags & CGLO_STOPCM;
302     bGStat   = flags & CGLO_GSTAT;
303
304     bReadEkin     = (flags & CGLO_READEKIN);
305     bScaleEkin    = (flags & CGLO_SCALEEKIN);
306     bEner         = flags & CGLO_ENERGY;
307     bTemp         = flags & CGLO_TEMPERATURE;
308     bPres         = (flags & CGLO_PRESSURE);
309     bConstrain    = (flags & CGLO_CONSTRAINT);
310     bIterate      = (flags & CGLO_ITERATE);
311     bFirstIterate = (flags & CGLO_FIRSTITERATE);
312
313     /* we calculate a full state kinetic energy either with full-step velocity verlet
314        or half step where we need the pressure */
315
316     bEkinAveVel = (ir->eI == eiVV || (ir->eI == eiVVAK && bPres) || bReadEkin);
317
318     /* in initalization, it sums the shake virial in vv, and to
319        sums ekinh_old in leapfrog (or if we are calculating ekinh_old) for other reasons */
320
321     /* ########## Kinetic energy  ############## */
322
323     if (bTemp)
324     {
325         /* Non-equilibrium MD: this is parallellized, but only does communication
326          * when there really is NEMD.
327          */
328
329         if (PAR(cr) && (ekind->bNEMD))
330         {
331             accumulate_u(cr, &(ir->opts), ekind);
332         }
333         debug_gmx();
334         if (bReadEkin)
335         {
336             restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
337         }
338         else
339         {
340
341             calc_ke_part(state, &(ir->opts), mdatoms, ekind, nrnb, bEkinAveVel, bIterate);
342         }
343
344         debug_gmx();
345     }
346
347     /* Calculate center of mass velocity if necessary, also parallellized */
348     if (bStopCM)
349     {
350         calc_vcm_grp(fplog, mdatoms->start, mdatoms->homenr, mdatoms,
351                      state->x, state->v, vcm);
352     }
353
354     if (bTemp || bStopCM || bPres || bEner || bConstrain)
355     {
356         if (!bGStat)
357         {
358             /* We will not sum ekinh_old,
359              * so signal that we still have to do it.
360              */
361             *bSumEkinhOld = TRUE;
362
363         }
364         else
365         {
366             if (gs != NULL)
367             {
368                 for (i = 0; i < eglsNR; i++)
369                 {
370                     gs_buf[i] = gs->sig[i];
371                 }
372             }
373             if (PAR(cr))
374             {
375                 wallcycle_start(wcycle, ewcMoveE);
376                 GMX_MPE_LOG(ev_global_stat_start);
377                 global_stat(fplog, gstat, cr, enerd, force_vir, shake_vir, mu_tot,
378                             ir, ekind, constr, bStopCM ? vcm : NULL,
379                             gs != NULL ? eglsNR : 0, gs_buf,
380                             top_global, state,
381                             *bSumEkinhOld, flags);
382                 GMX_MPE_LOG(ev_global_stat_finish);
383                 wallcycle_stop(wcycle, ewcMoveE);
384             }
385             if (gs != NULL)
386             {
387                 if (MULTISIM(cr) && bInterSimGS)
388                 {
389                     if (MASTER(cr))
390                     {
391                         /* Communicate the signals between the simulations */
392                         gmx_sum_sim(eglsNR, gs_buf, cr->ms);
393                     }
394                     /* Communicate the signals form the master to the others */
395                     gmx_bcast(eglsNR*sizeof(gs_buf[0]), gs_buf, cr);
396                 }
397                 for (i = 0; i < eglsNR; i++)
398                 {
399                     if (bInterSimGS || gs_simlocal[i])
400                     {
401                         /* Set the communicated signal only when it is non-zero,
402                          * since signals might not be processed at each MD step.
403                          */
404                         gsi = (gs_buf[i] >= 0 ?
405                                (int)(gs_buf[i] + 0.5) :
406                                (int)(gs_buf[i] - 0.5));
407                         if (gsi != 0)
408                         {
409                             gs->set[i] = gsi;
410                         }
411                         /* Turn off the local signal */
412                         gs->sig[i] = 0;
413                     }
414                 }
415             }
416             *bSumEkinhOld = FALSE;
417         }
418     }
419
420     if (!ekind->bNEMD && debug && bTemp && (vcm->nr > 0))
421     {
422         correct_ekin(debug,
423                      mdatoms->start, mdatoms->start+mdatoms->homenr,
424                      state->v, vcm->group_p[0],
425                      mdatoms->massT, mdatoms->tmass, ekind->ekin);
426     }
427
428     /* Do center of mass motion removal */
429     if (bStopCM)
430     {
431         check_cm_grp(fplog, vcm, ir, 1);
432         do_stopcm_grp(fplog, mdatoms->start, mdatoms->homenr, mdatoms->cVCM,
433                       state->x, state->v, vcm);
434         inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
435     }
436
437     if (bEner)
438     {
439         /* Calculate the amplitude of the cosine velocity profile */
440         ekind->cosacc.vcos = ekind->cosacc.mvcos/mdatoms->tmass;
441     }
442
443     if (bTemp)
444     {
445         /* Sum the kinetic energies of the groups & calc temp */
446         /* compute full step kinetic energies if vv, or if vv-avek and we are computing the pressure with IR_NPT_TROTTER */
447         /* three maincase:  VV with AveVel (md-vv), vv with AveEkin (md-vv-avek), leap with AveEkin (md).
448            Leap with AveVel is not supported; it's not clear that it will actually work.
449            bEkinAveVel: If TRUE, we simply multiply ekin by ekinscale to get a full step kinetic energy.
450            If FALSE, we average ekinh_old and ekinh*ekinscale_nhc to get an averaged half step kinetic energy.
451            bSaveEkinOld: If TRUE (in the case of iteration = bIterate is TRUE), we don't reset the ekinscale_nhc.
452            If FALSE, we go ahead and erase over it.
453          */
454         enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, &dvdl_ekin,
455                                        bEkinAveVel, bIterate, bScaleEkin);
456         enerd->dvdl_lin[efptMASS] = (double) dvdl_ekin;
457
458         enerd->term[F_EKIN] = trace(ekind->ekin);
459     }
460
461     /* ##########  Long range energy information ###### */
462
463     if (bEner || bPres || bConstrain)
464     {
465         calc_dispcorr(fplog, ir, fr, 0, top_global->natoms, box, state->lambda[efptVDW],
466                       corr_pres, corr_vir, &prescorr, &enercorr, &dvdlcorr);
467     }
468
469     if (bEner && bFirstIterate)
470     {
471         enerd->term[F_DISPCORR]  = enercorr;
472         enerd->term[F_EPOT]     += enercorr;
473         enerd->term[F_DVDL_VDW] += dvdlcorr;
474     }
475
476     /* ########## Now pressure ############## */
477     if (bPres || bConstrain)
478     {
479
480         m_add(force_vir, shake_vir, total_vir);
481
482         /* Calculate pressure and apply LR correction if PPPM is used.
483          * Use the box from last timestep since we already called update().
484          */
485
486         enerd->term[F_PRES] = calc_pres(fr->ePBC, ir->nwall, box, ekind->ekin, total_vir, pres);
487
488         /* Calculate long range corrections to pressure and energy */
489         /* this adds to enerd->term[F_PRES] and enerd->term[F_ETOT],
490            and computes enerd->term[F_DISPCORR].  Also modifies the
491            total_vir and pres tesors */
492
493         m_add(total_vir, corr_vir, total_vir);
494         m_add(pres, corr_pres, pres);
495         enerd->term[F_PDISPCORR] = prescorr;
496         enerd->term[F_PRES]     += prescorr;
497         *pcurr                   = enerd->term[F_PRES];
498     }
499 }
500
501 void check_nst_param(FILE *fplog, t_commrec *cr,
502                      const char *desc_nst, int nst,
503                      const char *desc_p, int *p)
504 {
505     if (*p > 0 && *p % nst != 0)
506     {
507         /* Round up to the next multiple of nst */
508         *p = ((*p)/nst + 1)*nst;
509         md_print_warn(cr, fplog,
510                       "NOTE: %s changes %s to %d\n", desc_nst, desc_p, *p);
511     }
512 }
513
514 void set_current_lambdas(gmx_large_int_t step, t_lambda *fepvals, gmx_bool bRerunMD,
515                          t_trxframe *rerun_fr, t_state *state_global, t_state *state, double lam0[])
516 /* find the current lambdas.  If rerunning, we either read in a state, or a lambda value,
517    requiring different logic. */
518 {
519     real frac;
520     int  i, fep_state = 0;
521     if (bRerunMD)
522     {
523         if (rerun_fr->bLambda)
524         {
525             if (fepvals->delta_lambda==0)
526             {
527                 state_global->lambda[efptFEP] = rerun_fr->lambda;
528                 for (i = 0; i < efptNR; i++)
529                 {
530                     if (i != efptFEP)
531                     {
532                         state->lambda[i] = state_global->lambda[i];
533                     }
534                 }
535             }
536             else
537             {
538                 /* find out between which two value of lambda we should be */
539                 frac      = (step*fepvals->delta_lambda);
540                 fep_state = floor(frac*fepvals->n_lambda);
541                 /* interpolate between this state and the next */
542                 /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
543                 frac = (frac*fepvals->n_lambda)-fep_state;
544                 for (i = 0; i < efptNR; i++)
545                 {
546                     state_global->lambda[i] = lam0[i] + (fepvals->all_lambda[i][fep_state]) +
547                         frac*(fepvals->all_lambda[i][fep_state+1]-fepvals->all_lambda[i][fep_state]);
548                 }
549             }
550         }
551         else if (rerun_fr->bFepState)
552         {
553             state_global->fep_state = rerun_fr->fep_state;
554             for (i = 0; i < efptNR; i++)
555             {
556                 state_global->lambda[i] = fepvals->all_lambda[i][fep_state];
557             }
558         }
559     }
560     else
561     {
562         if (fepvals->delta_lambda != 0)
563         {
564             /* find out between which two value of lambda we should be */
565             frac = (step*fepvals->delta_lambda);
566             if (fepvals->n_lambda > 0)
567             {
568                 fep_state = floor(frac*fepvals->n_lambda);
569                 /* interpolate between this state and the next */
570                 /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
571                 frac = (frac*fepvals->n_lambda)-fep_state;
572                 for (i = 0; i < efptNR; i++)
573                 {
574                     state_global->lambda[i] = lam0[i] + (fepvals->all_lambda[i][fep_state]) +
575                         frac*(fepvals->all_lambda[i][fep_state+1]-fepvals->all_lambda[i][fep_state]);
576                 }
577             }
578             else
579             {
580                 for (i = 0; i < efptNR; i++)
581                 {
582                     state_global->lambda[i] = lam0[i] + frac;
583                 }
584             }
585         }
586         else
587         {
588             if (state->fep_state > 0) {
589                 state_global->fep_state = state->fep_state; /* state->fep is the one updated by bExpanded */
590                 for (i = 0; i < efptNR; i++)
591                 {
592                     state_global->lambda[i] = fepvals->all_lambda[i][state_global->fep_state];
593                 }
594             }
595         }
596     }
597     for (i = 0; i < efptNR; i++)
598     {
599         state->lambda[i] = state_global->lambda[i];
600     }
601 }
602
603 static void min_zero(int *n, int i)
604 {
605     if (i > 0 && (*n == 0 || i < *n))
606     {
607         *n = i;
608     }
609 }
610
611 static int lcd4(int i1, int i2, int i3, int i4)
612 {
613     int nst;
614
615     nst = 0;
616     min_zero(&nst, i1);
617     min_zero(&nst, i2);
618     min_zero(&nst, i3);
619     min_zero(&nst, i4);
620     if (nst == 0)
621     {
622         gmx_incons("All 4 inputs for determininig nstglobalcomm are <= 0");
623     }
624
625     while (nst > 1 && ((i1 > 0 && i1 % nst != 0)  ||
626                        (i2 > 0 && i2 % nst != 0)  ||
627                        (i3 > 0 && i3 % nst != 0)  ||
628                        (i4 > 0 && i4 % nst != 0)))
629     {
630         nst--;
631     }
632
633     return nst;
634 }
635
636 int check_nstglobalcomm(FILE *fplog, t_commrec *cr,
637                         int nstglobalcomm, t_inputrec *ir)
638 {
639     if (!EI_DYNAMICS(ir->eI))
640     {
641         nstglobalcomm = 1;
642     }
643
644     if (nstglobalcomm == -1)
645     {
646         if (!(ir->nstcalcenergy > 0 ||
647               ir->nstlist > 0 ||
648               ir->etc != etcNO ||
649               ir->epc != epcNO))
650         {
651             nstglobalcomm = 10;
652             if (ir->nstenergy > 0 && ir->nstenergy < nstglobalcomm)
653             {
654                 nstglobalcomm = ir->nstenergy;
655             }
656         }
657         else
658         {
659             /* Ensure that we do timely global communication for
660              * (possibly) each of the four following options.
661              */
662             nstglobalcomm = lcd4(ir->nstcalcenergy,
663                                  ir->nstlist,
664                                  ir->etc != etcNO ? ir->nsttcouple : 0,
665                                  ir->epc != epcNO ? ir->nstpcouple : 0);
666         }
667     }
668     else
669     {
670         if (ir->nstlist > 0 &&
671             nstglobalcomm > ir->nstlist && nstglobalcomm % ir->nstlist != 0)
672         {
673             nstglobalcomm = (nstglobalcomm / ir->nstlist)*ir->nstlist;
674             md_print_warn(cr, fplog, "WARNING: nstglobalcomm is larger than nstlist, but not a multiple, setting it to %d\n", nstglobalcomm);
675         }
676         if (ir->nstcalcenergy > 0)
677         {
678             check_nst_param(fplog, cr, "-gcom", nstglobalcomm,
679                             "nstcalcenergy", &ir->nstcalcenergy);
680         }
681         if (ir->etc != etcNO && ir->nsttcouple > 0)
682         {
683             check_nst_param(fplog, cr, "-gcom", nstglobalcomm,
684                             "nsttcouple", &ir->nsttcouple);
685         }
686         if (ir->epc != epcNO && ir->nstpcouple > 0)
687         {
688             check_nst_param(fplog, cr, "-gcom", nstglobalcomm,
689                             "nstpcouple", &ir->nstpcouple);
690         }
691
692         check_nst_param(fplog, cr, "-gcom", nstglobalcomm,
693                         "nstenergy", &ir->nstenergy);
694
695         check_nst_param(fplog, cr, "-gcom", nstglobalcomm,
696                         "nstlog", &ir->nstlog);
697     }
698
699     if (ir->comm_mode != ecmNO && ir->nstcomm < nstglobalcomm)
700     {
701         md_print_warn(cr, fplog, "WARNING: Changing nstcomm from %d to %d\n",
702                       ir->nstcomm, nstglobalcomm);
703         ir->nstcomm = nstglobalcomm;
704     }
705
706     return nstglobalcomm;
707 }
708
709 void check_ir_old_tpx_versions(t_commrec *cr, FILE *fplog,
710                                t_inputrec *ir, gmx_mtop_t *mtop)
711 {
712     /* Check required for old tpx files */
713     if (IR_TWINRANGE(*ir) && ir->nstlist > 1 &&
714         ir->nstcalcenergy % ir->nstlist != 0)
715     {
716         md_print_warn(cr, fplog, "Old tpr file with twin-range settings: modifying energy calculation and/or T/P-coupling frequencies\n");
717
718         if (gmx_mtop_ftype_count(mtop, F_CONSTR) +
719             gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0 &&
720             ir->eConstrAlg == econtSHAKE)
721         {
722             md_print_warn(cr, fplog, "With twin-range cut-off's and SHAKE the virial and pressure are incorrect\n");
723             if (ir->epc != epcNO)
724             {
725                 gmx_fatal(FARGS, "Can not do pressure coupling with twin-range cut-off's and SHAKE");
726             }
727         }
728         check_nst_param(fplog, cr, "nstlist", ir->nstlist,
729                         "nstcalcenergy", &ir->nstcalcenergy);
730         if (ir->epc != epcNO)
731         {
732             check_nst_param(fplog, cr, "nstlist", ir->nstlist,
733                             "nstpcouple", &ir->nstpcouple);
734         }
735         check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
736                         "nstenergy", &ir->nstenergy);
737         check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
738                         "nstlog", &ir->nstlog);
739         if (ir->efep != efepNO)
740         {
741             check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
742                             "nstdhdl", &ir->fepvals->nstdhdl);
743         }
744     }
745
746     if (EI_VV(ir->eI) && IR_TWINRANGE(*ir) && ir->nstlist > 1)
747     {
748         gmx_fatal(FARGS, "Twin-range multiple time stepping does not work with integrator %s.", ei_names[ir->eI]);
749     }
750 }
751
752 void rerun_parallel_comm(t_commrec *cr, t_trxframe *fr,
753                          gmx_bool *bNotLastFrame)
754 {
755     gmx_bool bAlloc;
756     rvec    *xp, *vp;
757
758     bAlloc = (fr->natoms == 0);
759
760     if (MASTER(cr) && !*bNotLastFrame)
761     {
762         fr->natoms = -1;
763     }
764     xp = fr->x;
765     vp = fr->v;
766     gmx_bcast(sizeof(*fr), fr, cr);
767     fr->x = xp;
768     fr->v = vp;
769
770     *bNotLastFrame = (fr->natoms >= 0);
771
772     if (*bNotLastFrame && PARTDECOMP(cr))
773     {
774         /* x and v are the only variable size quantities stored in trr
775          * that are required for rerun (f is not needed).
776          */
777         if (bAlloc)
778         {
779             snew(fr->x, fr->natoms);
780             snew(fr->v, fr->natoms);
781         }
782         if (fr->bX)
783         {
784             gmx_bcast(fr->natoms*sizeof(fr->x[0]), fr->x[0], cr);
785         }
786         if (fr->bV)
787         {
788             gmx_bcast(fr->natoms*sizeof(fr->v[0]), fr->v[0], cr);
789         }
790     }
791 }