Reduced usage of typdefs.h
[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 "mdrun.h"
43 #include "domdec.h"
44 #include "mtop_util.h"
45 #include "vcm.h"
46 #include "nrnb.h"
47 #include "macros.h"
48 #include "md_logging.h"
49 #include "md_support.h"
50
51 #include "gromacs/legacyheaders/types/commrec.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/timing/wallcycle.h"
54 #include "gromacs/utility/cstringutil.h"
55 #include "gromacs/utility/smalloc.h"
56
57 /* Is the signal in one simulation independent of other simulations? */
58 gmx_bool gs_simlocal[eglsNR] = { TRUE, FALSE, FALSE, TRUE };
59
60 /* check which of the multisim simulations has the shortest number of
61    steps and return that number of nsteps */
62 gmx_int64_t get_multisim_nsteps(const t_commrec *cr,
63                                 gmx_int64_t      nsteps)
64 {
65     gmx_int64_t steps_out;
66
67     if (MASTER(cr))
68     {
69         gmx_int64_t     *buf;
70         int              s;
71
72         snew(buf, cr->ms->nsim);
73
74         buf[cr->ms->sim] = nsteps;
75         gmx_sumli_sim(cr->ms->nsim, buf, cr->ms);
76
77         steps_out = -1;
78         for (s = 0; s < cr->ms->nsim; s++)
79         {
80             /* find the smallest positive number */
81             if (buf[s] >= 0 && ((steps_out < 0) || (buf[s] < steps_out)) )
82             {
83                 steps_out = buf[s];
84             }
85         }
86         sfree(buf);
87
88         /* if we're the limiting simulation, don't do anything */
89         if (steps_out >= 0 && steps_out < nsteps)
90         {
91             char strbuf[255];
92             snprintf(strbuf, 255, "Will stop simulation %%d after %s steps (another simulation will end then).\n", "%"GMX_PRId64);
93             fprintf(stderr, strbuf, cr->ms->sim, steps_out);
94         }
95     }
96     /* broadcast to non-masters */
97     gmx_bcast(sizeof(gmx_int64_t), &steps_out, cr);
98     return steps_out;
99 }
100
101 int multisim_min(const gmx_multisim_t *ms, int nmin, int n)
102 {
103     int     *buf;
104     gmx_bool bPos, bEqual;
105     int      s, d;
106
107     snew(buf, ms->nsim);
108     buf[ms->sim] = n;
109     gmx_sumi_sim(ms->nsim, buf, ms);
110     bPos   = TRUE;
111     bEqual = TRUE;
112     for (s = 0; s < ms->nsim; s++)
113     {
114         bPos   = bPos   && (buf[s] > 0);
115         bEqual = bEqual && (buf[s] == buf[0]);
116     }
117     if (bPos)
118     {
119         if (bEqual)
120         {
121             nmin = min(nmin, buf[0]);
122         }
123         else
124         {
125             /* Find the least common multiple */
126             for (d = 2; d < nmin; d++)
127             {
128                 s = 0;
129                 while (s < ms->nsim && d % buf[s] == 0)
130                 {
131                     s++;
132                 }
133                 if (s == ms->nsim)
134                 {
135                     /* We found the LCM and it is less than nmin */
136                     nmin = d;
137                     break;
138                 }
139             }
140         }
141     }
142     sfree(buf);
143
144     return nmin;
145 }
146
147 int multisim_nstsimsync(const t_commrec *cr,
148                         const t_inputrec *ir, int repl_ex_nst)
149 {
150     int nmin;
151
152     if (MASTER(cr))
153     {
154         nmin = INT_MAX;
155         nmin = multisim_min(cr->ms, nmin, ir->nstlist);
156         nmin = multisim_min(cr->ms, nmin, ir->nstcalcenergy);
157         nmin = multisim_min(cr->ms, nmin, repl_ex_nst);
158         if (nmin == INT_MAX)
159         {
160             gmx_fatal(FARGS, "Can not find an appropriate interval for inter-simulation communication, since nstlist, nstcalcenergy and -replex are all <= 0");
161         }
162         /* Avoid inter-simulation communication at every (second) step */
163         if (nmin <= 2)
164         {
165             nmin = 10;
166         }
167     }
168
169     gmx_bcast(sizeof(int), &nmin, cr);
170
171     return nmin;
172 }
173
174 void init_global_signals(globsig_t *gs, const t_commrec *cr,
175                          const t_inputrec *ir, int repl_ex_nst)
176 {
177     int i;
178
179     if (MULTISIM(cr))
180     {
181         gs->nstms = multisim_nstsimsync(cr, ir, repl_ex_nst);
182         if (debug)
183         {
184             fprintf(debug, "Syncing simulations for checkpointing and termination every %d steps\n", gs->nstms);
185         }
186     }
187     else
188     {
189         gs->nstms = 1;
190     }
191
192     for (i = 0; i < eglsNR; i++)
193     {
194         gs->sig[i] = 0;
195         gs->set[i] = 0;
196     }
197 }
198
199 void copy_coupling_state(t_state *statea, t_state *stateb,
200                          gmx_ekindata_t *ekinda, gmx_ekindata_t *ekindb, t_grpopts* opts)
201 {
202
203     /* MRS note -- might be able to get rid of some of the arguments.  Look over it when it's all debugged */
204
205     int i, j, nc;
206
207     /* Make sure we have enough space for x and v */
208     if (statea->nalloc > stateb->nalloc)
209     {
210         stateb->nalloc = statea->nalloc;
211         srenew(stateb->x, stateb->nalloc);
212         srenew(stateb->v, stateb->nalloc);
213     }
214
215     stateb->natoms     = statea->natoms;
216     stateb->ngtc       = statea->ngtc;
217     stateb->nnhpres    = statea->nnhpres;
218     stateb->veta       = statea->veta;
219     if (ekinda)
220     {
221         copy_mat(ekinda->ekin, ekindb->ekin);
222         for (i = 0; i < stateb->ngtc; i++)
223         {
224             ekindb->tcstat[i].T  = ekinda->tcstat[i].T;
225             ekindb->tcstat[i].Th = ekinda->tcstat[i].Th;
226             copy_mat(ekinda->tcstat[i].ekinh, ekindb->tcstat[i].ekinh);
227             copy_mat(ekinda->tcstat[i].ekinf, ekindb->tcstat[i].ekinf);
228             ekindb->tcstat[i].ekinscalef_nhc =  ekinda->tcstat[i].ekinscalef_nhc;
229             ekindb->tcstat[i].ekinscaleh_nhc =  ekinda->tcstat[i].ekinscaleh_nhc;
230             ekindb->tcstat[i].vscale_nhc     =  ekinda->tcstat[i].vscale_nhc;
231         }
232     }
233     copy_rvecn(statea->x, stateb->x, 0, stateb->natoms);
234     copy_rvecn(statea->v, stateb->v, 0, stateb->natoms);
235     copy_mat(statea->box, stateb->box);
236     copy_mat(statea->box_rel, stateb->box_rel);
237     copy_mat(statea->boxv, stateb->boxv);
238
239     for (i = 0; i < stateb->ngtc; i++)
240     {
241         nc = i*opts->nhchainlength;
242         for (j = 0; j < opts->nhchainlength; j++)
243         {
244             stateb->nosehoover_xi[nc+j]  = statea->nosehoover_xi[nc+j];
245             stateb->nosehoover_vxi[nc+j] = statea->nosehoover_vxi[nc+j];
246         }
247     }
248     if (stateb->nhpres_xi != NULL)
249     {
250         for (i = 0; i < stateb->nnhpres; i++)
251         {
252             nc = i*opts->nhchainlength;
253             for (j = 0; j < opts->nhchainlength; j++)
254             {
255                 stateb->nhpres_xi[nc+j]  = statea->nhpres_xi[nc+j];
256                 stateb->nhpres_vxi[nc+j] = statea->nhpres_vxi[nc+j];
257             }
258         }
259     }
260 }
261
262 real compute_conserved_from_auxiliary(t_inputrec *ir, t_state *state, t_extmass *MassQ)
263 {
264     real quantity = 0;
265     switch (ir->etc)
266     {
267         case etcNO:
268             break;
269         case etcBERENDSEN:
270             break;
271         case etcNOSEHOOVER:
272             quantity = NPT_energy(ir, state, MassQ);
273             break;
274         case etcVRESCALE:
275             quantity = vrescale_energy(&(ir->opts), state->therm_integral);
276             break;
277         default:
278             break;
279     }
280     return quantity;
281 }
282
283 void compute_globals(FILE *fplog, gmx_global_stat_t gstat, t_commrec *cr, t_inputrec *ir,
284                      t_forcerec *fr, gmx_ekindata_t *ekind,
285                      t_state *state, t_state *state_global, t_mdatoms *mdatoms,
286                      t_nrnb *nrnb, t_vcm *vcm, gmx_wallcycle_t wcycle,
287                      gmx_enerdata_t *enerd, tensor force_vir, tensor shake_vir, tensor total_vir,
288                      tensor pres, rvec mu_tot, gmx_constr_t constr,
289                      globsig_t *gs, gmx_bool bInterSimGS,
290                      matrix box, gmx_mtop_t *top_global,
291                      gmx_bool *bSumEkinhOld, int flags)
292 {
293     int      i, gsi;
294     real     gs_buf[eglsNR];
295     tensor   corr_vir, corr_pres;
296     gmx_bool bEner, bPres, bTemp, bVV;
297     gmx_bool bRerunMD, bStopCM, bGStat, bIterate,
298              bFirstIterate, bReadEkin, bEkinAveVel, bScaleEkin, bConstrain;
299     real     ekin, temp, prescorr, enercorr, dvdlcorr, dvdl_ekin;
300
301     /* translate CGLO flags to gmx_booleans */
302     bRerunMD = flags & CGLO_RERUNMD;
303     bStopCM  = flags & CGLO_STOPCM;
304     bGStat   = flags & CGLO_GSTAT;
305
306     bReadEkin     = (flags & CGLO_READEKIN);
307     bScaleEkin    = (flags & CGLO_SCALEEKIN);
308     bEner         = flags & CGLO_ENERGY;
309     bTemp         = flags & CGLO_TEMPERATURE;
310     bPres         = (flags & CGLO_PRESSURE);
311     bConstrain    = (flags & CGLO_CONSTRAINT);
312     bIterate      = (flags & CGLO_ITERATE);
313     bFirstIterate = (flags & CGLO_FIRSTITERATE);
314
315     /* we calculate a full state kinetic energy either with full-step velocity verlet
316        or half step where we need the pressure */
317
318     bEkinAveVel = (ir->eI == eiVV || (ir->eI == eiVVAK && bPres) || bReadEkin);
319
320     /* in initalization, it sums the shake virial in vv, and to
321        sums ekinh_old in leapfrog (or if we are calculating ekinh_old) for other reasons */
322
323     /* ########## Kinetic energy  ############## */
324
325     if (bTemp)
326     {
327         /* Non-equilibrium MD: this is parallellized, but only does communication
328          * when there really is NEMD.
329          */
330
331         if (PAR(cr) && (ekind->bNEMD))
332         {
333             accumulate_u(cr, &(ir->opts), ekind);
334         }
335         debug_gmx();
336         if (bReadEkin)
337         {
338             restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
339         }
340         else
341         {
342
343             calc_ke_part(state, &(ir->opts), mdatoms, ekind, nrnb, bEkinAveVel, bIterate);
344         }
345
346         debug_gmx();
347     }
348
349     /* Calculate center of mass velocity if necessary, also parallellized */
350     if (bStopCM)
351     {
352         calc_vcm_grp(0, mdatoms->homenr, mdatoms,
353                      state->x, state->v, vcm);
354     }
355
356     if (bTemp || bStopCM || bPres || bEner || bConstrain)
357     {
358         if (!bGStat)
359         {
360             /* We will not sum ekinh_old,
361              * so signal that we still have to do it.
362              */
363             *bSumEkinhOld = TRUE;
364
365         }
366         else
367         {
368             if (gs != NULL)
369             {
370                 for (i = 0; i < eglsNR; i++)
371                 {
372                     gs_buf[i] = gs->sig[i];
373                 }
374             }
375             if (PAR(cr))
376             {
377                 wallcycle_start(wcycle, ewcMoveE);
378                 global_stat(fplog, gstat, cr, enerd, force_vir, shake_vir, mu_tot,
379                             ir, ekind, constr, bStopCM ? vcm : NULL,
380                             gs != NULL ? eglsNR : 0, gs_buf,
381                             top_global, state,
382                             *bSumEkinhOld, flags);
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                      0, 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(0, 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, 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     }
498 }
499
500 void check_nst_param(FILE *fplog, t_commrec *cr,
501                      const char *desc_nst, int nst,
502                      const char *desc_p, int *p)
503 {
504     if (*p > 0 && *p % nst != 0)
505     {
506         /* Round up to the next multiple of nst */
507         *p = ((*p)/nst + 1)*nst;
508         md_print_warn(cr, fplog,
509                       "NOTE: %s changes %s to %d\n", desc_nst, desc_p, *p);
510     }
511 }
512
513 void set_current_lambdas(gmx_int64_t step, t_lambda *fepvals, gmx_bool bRerunMD,
514                          t_trxframe *rerun_fr, t_state *state_global, t_state *state, double lam0[])
515 /* find the current lambdas.  If rerunning, we either read in a state, or a lambda value,
516    requiring different logic. */
517 {
518     real frac;
519     int  i, fep_state = 0;
520     if (bRerunMD)
521     {
522         if (rerun_fr->bLambda)
523         {
524             if (fepvals->delta_lambda == 0)
525             {
526                 state_global->lambda[efptFEP] = rerun_fr->lambda;
527                 for (i = 0; i < efptNR; i++)
528                 {
529                     if (i != efptFEP)
530                     {
531                         state->lambda[i] = state_global->lambda[i];
532                     }
533                 }
534             }
535             else
536             {
537                 /* find out between which two value of lambda we should be */
538                 frac      = (step*fepvals->delta_lambda);
539                 fep_state = floor(frac*fepvals->n_lambda);
540                 /* interpolate between this state and the next */
541                 /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
542                 frac = (frac*fepvals->n_lambda)-fep_state;
543                 for (i = 0; i < efptNR; i++)
544                 {
545                     state_global->lambda[i] = lam0[i] + (fepvals->all_lambda[i][fep_state]) +
546                         frac*(fepvals->all_lambda[i][fep_state+1]-fepvals->all_lambda[i][fep_state]);
547                 }
548             }
549         }
550         else if (rerun_fr->bFepState)
551         {
552             state_global->fep_state = rerun_fr->fep_state;
553             for (i = 0; i < efptNR; i++)
554             {
555                 state_global->lambda[i] = fepvals->all_lambda[i][fep_state];
556             }
557         }
558     }
559     else
560     {
561         if (fepvals->delta_lambda != 0)
562         {
563             /* find out between which two value of lambda we should be */
564             frac = (step*fepvals->delta_lambda);
565             if (fepvals->n_lambda > 0)
566             {
567                 fep_state = floor(frac*fepvals->n_lambda);
568                 /* interpolate between this state and the next */
569                 /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
570                 frac = (frac*fepvals->n_lambda)-fep_state;
571                 for (i = 0; i < efptNR; i++)
572                 {
573                     state_global->lambda[i] = lam0[i] + (fepvals->all_lambda[i][fep_state]) +
574                         frac*(fepvals->all_lambda[i][fep_state+1]-fepvals->all_lambda[i][fep_state]);
575                 }
576             }
577             else
578             {
579                 for (i = 0; i < efptNR; i++)
580                 {
581                     state_global->lambda[i] = lam0[i] + frac;
582                 }
583             }
584         }
585         else
586         {
587             if (state->fep_state > 0)
588             {
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
747 void rerun_parallel_comm(t_commrec *cr, t_trxframe *fr,
748                          gmx_bool *bNotLastFrame)
749 {
750     gmx_bool bAlloc;
751     rvec    *xp, *vp;
752
753     bAlloc = (fr->natoms == 0);
754
755     if (MASTER(cr) && !*bNotLastFrame)
756     {
757         fr->natoms = -1;
758     }
759     xp = fr->x;
760     vp = fr->v;
761     gmx_bcast(sizeof(*fr), fr, cr);
762     fr->x = xp;
763     fr->v = vp;
764
765     *bNotLastFrame = (fr->natoms >= 0);
766
767 }