Simplified section "Porting".
[alexxy/gromacs.git] / src / gromacs / 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  * Copyright (c) 2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include "typedefs.h"
42 #include "gromacs/utility/cstringutil.h"
43 #include "gromacs/utility/smalloc.h"
44 #include "mdrun.h"
45 #include "domdec.h"
46 #include "mtop_util.h"
47 #include "vcm.h"
48 #include "nrnb.h"
49 #include "macros.h"
50 #include "md_logging.h"
51 #include "md_support.h"
52
53 #include "gromacs/timing/wallcycle.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_int64_t get_multisim_nsteps(const t_commrec *cr,
61                                 gmx_int64_t      nsteps)
62 {
63     gmx_int64_t steps_out;
64
65     if (MASTER(cr))
66     {
67         gmx_int64_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_PRId64);
91             fprintf(stderr, strbuf, cr->ms->sim, steps_out);
92         }
93     }
94     /* broadcast to non-masters */
95     gmx_bcast(sizeof(gmx_int64_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,
289                      gmx_bool *bSumEkinhOld, int flags)
290 {
291     int      i, gsi;
292     real     gs_buf[eglsNR];
293     tensor   corr_vir, corr_pres;
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(0, 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                 global_stat(fplog, gstat, cr, enerd, force_vir, shake_vir, mu_tot,
377                             ir, ekind, constr, bStopCM ? vcm : NULL,
378                             gs != NULL ? eglsNR : 0, gs_buf,
379                             top_global, state,
380                             *bSumEkinhOld, flags);
381                 wallcycle_stop(wcycle, ewcMoveE);
382             }
383             if (gs != NULL)
384             {
385                 if (MULTISIM(cr) && bInterSimGS)
386                 {
387                     if (MASTER(cr))
388                     {
389                         /* Communicate the signals between the simulations */
390                         gmx_sum_sim(eglsNR, gs_buf, cr->ms);
391                     }
392                     /* Communicate the signals form the master to the others */
393                     gmx_bcast(eglsNR*sizeof(gs_buf[0]), gs_buf, cr);
394                 }
395                 for (i = 0; i < eglsNR; i++)
396                 {
397                     if (bInterSimGS || gs_simlocal[i])
398                     {
399                         /* Set the communicated signal only when it is non-zero,
400                          * since signals might not be processed at each MD step.
401                          */
402                         gsi = (gs_buf[i] >= 0 ?
403                                (int)(gs_buf[i] + 0.5) :
404                                (int)(gs_buf[i] - 0.5));
405                         if (gsi != 0)
406                         {
407                             gs->set[i] = gsi;
408                         }
409                         /* Turn off the local signal */
410                         gs->sig[i] = 0;
411                     }
412                 }
413             }
414             *bSumEkinhOld = FALSE;
415         }
416     }
417
418     if (!ekind->bNEMD && debug && bTemp && (vcm->nr > 0))
419     {
420         correct_ekin(debug,
421                      0, mdatoms->homenr,
422                      state->v, vcm->group_p[0],
423                      mdatoms->massT, mdatoms->tmass, ekind->ekin);
424     }
425
426     /* Do center of mass motion removal */
427     if (bStopCM)
428     {
429         check_cm_grp(fplog, vcm, ir, 1);
430         do_stopcm_grp(0, mdatoms->homenr, mdatoms->cVCM,
431                       state->x, state->v, vcm);
432         inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
433     }
434
435     if (bEner)
436     {
437         /* Calculate the amplitude of the cosine velocity profile */
438         ekind->cosacc.vcos = ekind->cosacc.mvcos/mdatoms->tmass;
439     }
440
441     if (bTemp)
442     {
443         /* Sum the kinetic energies of the groups & calc temp */
444         /* compute full step kinetic energies if vv, or if vv-avek and we are computing the pressure with IR_NPT_TROTTER */
445         /* three maincase:  VV with AveVel (md-vv), vv with AveEkin (md-vv-avek), leap with AveEkin (md).
446            Leap with AveVel is not supported; it's not clear that it will actually work.
447            bEkinAveVel: If TRUE, we simply multiply ekin by ekinscale to get a full step kinetic energy.
448            If FALSE, we average ekinh_old and ekinh*ekinscale_nhc to get an averaged half step kinetic energy.
449            bSaveEkinOld: If TRUE (in the case of iteration = bIterate is TRUE), we don't reset the ekinscale_nhc.
450            If FALSE, we go ahead and erase over it.
451          */
452         enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, &dvdl_ekin,
453                                        bEkinAveVel, bScaleEkin);
454         enerd->dvdl_lin[efptMASS] = (double) dvdl_ekin;
455
456         enerd->term[F_EKIN] = trace(ekind->ekin);
457     }
458
459     /* ##########  Long range energy information ###### */
460
461     if (bEner || bPres || bConstrain)
462     {
463         calc_dispcorr(fplog, ir, fr, 0, top_global->natoms, box, state->lambda[efptVDW],
464                       corr_pres, corr_vir, &prescorr, &enercorr, &dvdlcorr);
465     }
466
467     if (bEner && bFirstIterate)
468     {
469         enerd->term[F_DISPCORR]  = enercorr;
470         enerd->term[F_EPOT]     += enercorr;
471         enerd->term[F_DVDL_VDW] += dvdlcorr;
472     }
473
474     /* ########## Now pressure ############## */
475     if (bPres || bConstrain)
476     {
477
478         m_add(force_vir, shake_vir, total_vir);
479
480         /* Calculate pressure and apply LR correction if PPPM is used.
481          * Use the box from last timestep since we already called update().
482          */
483
484         enerd->term[F_PRES] = calc_pres(fr->ePBC, ir->nwall, box, ekind->ekin, total_vir, pres);
485
486         /* Calculate long range corrections to pressure and energy */
487         /* this adds to enerd->term[F_PRES] and enerd->term[F_ETOT],
488            and computes enerd->term[F_DISPCORR].  Also modifies the
489            total_vir and pres tesors */
490
491         m_add(total_vir, corr_vir, total_vir);
492         m_add(pres, corr_pres, pres);
493         enerd->term[F_PDISPCORR] = prescorr;
494         enerd->term[F_PRES]     += prescorr;
495     }
496 }
497
498 void check_nst_param(FILE *fplog, t_commrec *cr,
499                      const char *desc_nst, int nst,
500                      const char *desc_p, int *p)
501 {
502     if (*p > 0 && *p % nst != 0)
503     {
504         /* Round up to the next multiple of nst */
505         *p = ((*p)/nst + 1)*nst;
506         md_print_warn(cr, fplog,
507                       "NOTE: %s changes %s to %d\n", desc_nst, desc_p, *p);
508     }
509 }
510
511 void set_current_lambdas(gmx_int64_t step, t_lambda *fepvals, gmx_bool bRerunMD,
512                          t_trxframe *rerun_fr, t_state *state_global, t_state *state, double lam0[])
513 /* find the current lambdas.  If rerunning, we either read in a state, or a lambda value,
514    requiring different logic. */
515 {
516     real frac;
517     int  i, fep_state = 0;
518     if (bRerunMD)
519     {
520         if (rerun_fr->bLambda)
521         {
522             if (fepvals->delta_lambda == 0)
523             {
524                 state_global->lambda[efptFEP] = rerun_fr->lambda;
525                 for (i = 0; i < efptNR; i++)
526                 {
527                     if (i != efptFEP)
528                     {
529                         state->lambda[i] = state_global->lambda[i];
530                     }
531                 }
532             }
533             else
534             {
535                 /* find out between which two value of lambda we should be */
536                 frac      = (step*fepvals->delta_lambda);
537                 fep_state = floor(frac*fepvals->n_lambda);
538                 /* interpolate between this state and the next */
539                 /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
540                 frac = (frac*fepvals->n_lambda)-fep_state;
541                 for (i = 0; i < efptNR; i++)
542                 {
543                     state_global->lambda[i] = lam0[i] + (fepvals->all_lambda[i][fep_state]) +
544                         frac*(fepvals->all_lambda[i][fep_state+1]-fepvals->all_lambda[i][fep_state]);
545                 }
546             }
547         }
548         else if (rerun_fr->bFepState)
549         {
550             state_global->fep_state = rerun_fr->fep_state;
551             for (i = 0; i < efptNR; i++)
552             {
553                 state_global->lambda[i] = fepvals->all_lambda[i][fep_state];
554             }
555         }
556     }
557     else
558     {
559         if (fepvals->delta_lambda != 0)
560         {
561             /* find out between which two value of lambda we should be */
562             frac = (step*fepvals->delta_lambda);
563             if (fepvals->n_lambda > 0)
564             {
565                 fep_state = floor(frac*fepvals->n_lambda);
566                 /* interpolate between this state and the next */
567                 /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
568                 frac = (frac*fepvals->n_lambda)-fep_state;
569                 for (i = 0; i < efptNR; i++)
570                 {
571                     state_global->lambda[i] = lam0[i] + (fepvals->all_lambda[i][fep_state]) +
572                         frac*(fepvals->all_lambda[i][fep_state+1]-fepvals->all_lambda[i][fep_state]);
573                 }
574             }
575             else
576             {
577                 for (i = 0; i < efptNR; i++)
578                 {
579                     state_global->lambda[i] = lam0[i] + frac;
580                 }
581             }
582         }
583         else
584         {
585             if (state->fep_state > 0)
586             {
587                 state_global->fep_state = state->fep_state; /* state->fep is the one updated by bExpanded */
588                 for (i = 0; i < efptNR; i++)
589                 {
590                     state_global->lambda[i] = fepvals->all_lambda[i][state_global->fep_state];
591                 }
592             }
593         }
594     }
595     for (i = 0; i < efptNR; i++)
596     {
597         state->lambda[i] = state_global->lambda[i];
598     }
599 }
600
601 static void min_zero(int *n, int i)
602 {
603     if (i > 0 && (*n == 0 || i < *n))
604     {
605         *n = i;
606     }
607 }
608
609 static int lcd4(int i1, int i2, int i3, int i4)
610 {
611     int nst;
612
613     nst = 0;
614     min_zero(&nst, i1);
615     min_zero(&nst, i2);
616     min_zero(&nst, i3);
617     min_zero(&nst, i4);
618     if (nst == 0)
619     {
620         gmx_incons("All 4 inputs for determininig nstglobalcomm are <= 0");
621     }
622
623     while (nst > 1 && ((i1 > 0 && i1 % nst != 0)  ||
624                        (i2 > 0 && i2 % nst != 0)  ||
625                        (i3 > 0 && i3 % nst != 0)  ||
626                        (i4 > 0 && i4 % nst != 0)))
627     {
628         nst--;
629     }
630
631     return nst;
632 }
633
634 int check_nstglobalcomm(FILE *fplog, t_commrec *cr,
635                         int nstglobalcomm, t_inputrec *ir)
636 {
637     if (!EI_DYNAMICS(ir->eI))
638     {
639         nstglobalcomm = 1;
640     }
641
642     if (nstglobalcomm == -1)
643     {
644         if (!(ir->nstcalcenergy > 0 ||
645               ir->nstlist > 0 ||
646               ir->etc != etcNO ||
647               ir->epc != epcNO))
648         {
649             nstglobalcomm = 10;
650             if (ir->nstenergy > 0 && ir->nstenergy < nstglobalcomm)
651             {
652                 nstglobalcomm = ir->nstenergy;
653             }
654         }
655         else
656         {
657             /* Ensure that we do timely global communication for
658              * (possibly) each of the four following options.
659              */
660             nstglobalcomm = lcd4(ir->nstcalcenergy,
661                                  ir->nstlist,
662                                  ir->etc != etcNO ? ir->nsttcouple : 0,
663                                  ir->epc != epcNO ? ir->nstpcouple : 0);
664         }
665     }
666     else
667     {
668         if (ir->nstlist > 0 &&
669             nstglobalcomm > ir->nstlist && nstglobalcomm % ir->nstlist != 0)
670         {
671             nstglobalcomm = (nstglobalcomm / ir->nstlist)*ir->nstlist;
672             md_print_warn(cr, fplog, "WARNING: nstglobalcomm is larger than nstlist, but not a multiple, setting it to %d\n", nstglobalcomm);
673         }
674         if (ir->nstcalcenergy > 0)
675         {
676             check_nst_param(fplog, cr, "-gcom", nstglobalcomm,
677                             "nstcalcenergy", &ir->nstcalcenergy);
678         }
679         if (ir->etc != etcNO && ir->nsttcouple > 0)
680         {
681             check_nst_param(fplog, cr, "-gcom", nstglobalcomm,
682                             "nsttcouple", &ir->nsttcouple);
683         }
684         if (ir->epc != epcNO && ir->nstpcouple > 0)
685         {
686             check_nst_param(fplog, cr, "-gcom", nstglobalcomm,
687                             "nstpcouple", &ir->nstpcouple);
688         }
689
690         check_nst_param(fplog, cr, "-gcom", nstglobalcomm,
691                         "nstenergy", &ir->nstenergy);
692
693         check_nst_param(fplog, cr, "-gcom", nstglobalcomm,
694                         "nstlog", &ir->nstlog);
695     }
696
697     if (ir->comm_mode != ecmNO && ir->nstcomm < nstglobalcomm)
698     {
699         md_print_warn(cr, fplog, "WARNING: Changing nstcomm from %d to %d\n",
700                       ir->nstcomm, nstglobalcomm);
701         ir->nstcomm = nstglobalcomm;
702     }
703
704     return nstglobalcomm;
705 }
706
707 void check_ir_old_tpx_versions(t_commrec *cr, FILE *fplog,
708                                t_inputrec *ir, gmx_mtop_t *mtop)
709 {
710     /* Check required for old tpx files */
711     if (IR_TWINRANGE(*ir) && ir->nstlist > 1 &&
712         ir->nstcalcenergy % ir->nstlist != 0)
713     {
714         md_print_warn(cr, fplog, "Old tpr file with twin-range settings: modifying energy calculation and/or T/P-coupling frequencies\n");
715
716         if (gmx_mtop_ftype_count(mtop, F_CONSTR) +
717             gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0 &&
718             ir->eConstrAlg == econtSHAKE)
719         {
720             md_print_warn(cr, fplog, "With twin-range cut-off's and SHAKE the virial and pressure are incorrect\n");
721             if (ir->epc != epcNO)
722             {
723                 gmx_fatal(FARGS, "Can not do pressure coupling with twin-range cut-off's and SHAKE");
724             }
725         }
726         check_nst_param(fplog, cr, "nstlist", ir->nstlist,
727                         "nstcalcenergy", &ir->nstcalcenergy);
728         if (ir->epc != epcNO)
729         {
730             check_nst_param(fplog, cr, "nstlist", ir->nstlist,
731                             "nstpcouple", &ir->nstpcouple);
732         }
733         check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
734                         "nstenergy", &ir->nstenergy);
735         check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
736                         "nstlog", &ir->nstlog);
737         if (ir->efep != efepNO)
738         {
739             check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
740                             "nstdhdl", &ir->fepvals->nstdhdl);
741         }
742     }
743 }
744
745 void rerun_parallel_comm(t_commrec *cr, t_trxframe *fr,
746                          gmx_bool *bNotLastFrame)
747 {
748     gmx_bool bAlloc;
749     rvec    *xp, *vp;
750
751     bAlloc = (fr->natoms == 0);
752
753     if (MASTER(cr) && !*bNotLastFrame)
754     {
755         fr->natoms = -1;
756     }
757     xp = fr->x;
758     vp = fr->v;
759     gmx_bcast(sizeof(*fr), fr, cr);
760     fr->x = xp;
761     fr->v = vp;
762
763     *bNotLastFrame = (fr->natoms >= 0);
764
765 }