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