Divide default communicator from DD communicators
[alexxy/gromacs.git] / src / gromacs / mdlib / md_support.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017 The GROMACS development team.
7  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source 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
39 #include "gmxpre.h"
40
41 #include "md_support.h"
42
43 #include <climits>
44 #include <cmath>
45
46 #include <algorithm>
47
48 #include "gromacs/domdec/domdec.h"
49 #include "gromacs/gmxlib/network.h"
50 #include "gromacs/gmxlib/nrnb.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/mdlib/coupling.h"
53 #include "gromacs/mdlib/gmx_omp_nthreads.h"
54 #include "gromacs/mdlib/simulationsignal.h"
55 #include "gromacs/mdlib/stat.h"
56 #include "gromacs/mdlib/tgroup.h"
57 #include "gromacs/mdlib/update.h"
58 #include "gromacs/mdlib/vcm.h"
59 #include "gromacs/mdrunutility/multisim.h"
60 #include "gromacs/mdtypes/commrec.h"
61 #include "gromacs/mdtypes/df_history.h"
62 #include "gromacs/mdtypes/enerdata.h"
63 #include "gromacs/mdtypes/energyhistory.h"
64 #include "gromacs/mdtypes/forcerec.h"
65 #include "gromacs/mdtypes/group.h"
66 #include "gromacs/mdtypes/inputrec.h"
67 #include "gromacs/mdtypes/md_enums.h"
68 #include "gromacs/mdtypes/mdatom.h"
69 #include "gromacs/mdtypes/state.h"
70 #include "gromacs/pbcutil/pbc.h"
71 #include "gromacs/pulling/pull.h"
72 #include "gromacs/timing/wallcycle.h"
73 #include "gromacs/topology/mtop_util.h"
74 #include "gromacs/trajectory/trajectoryframe.h"
75 #include "gromacs/utility/arrayref.h"
76 #include "gromacs/utility/cstringutil.h"
77 #include "gromacs/utility/fatalerror.h"
78 #include "gromacs/utility/gmxassert.h"
79 #include "gromacs/utility/logger.h"
80 #include "gromacs/utility/smalloc.h"
81 #include "gromacs/utility/snprintf.h"
82
83 // TODO move this to multi-sim module
84 bool multisim_int_all_are_equal(const gmx_multisim_t* ms, int64_t value)
85 {
86     bool     allValuesAreEqual = true;
87     int64_t* buf;
88
89     GMX_RELEASE_ASSERT(ms, "Invalid use of multi-simulation pointer");
90
91     snew(buf, ms->numSimulations_);
92     /* send our value to all other master ranks, receive all of theirs */
93     buf[ms->simulationIndex_] = value;
94     gmx_sumli_sim(ms->numSimulations_, buf, ms);
95
96     for (int s = 0; s < ms->numSimulations_; s++)
97     {
98         if (buf[s] != value)
99         {
100             allValuesAreEqual = false;
101             break;
102         }
103     }
104
105     sfree(buf);
106
107     return allValuesAreEqual;
108 }
109
110 int multisim_min(const gmx_multisim_t* ms, int nmin, int n)
111 {
112     int*     buf;
113     gmx_bool bPos, bEqual;
114     int      s, d;
115
116     snew(buf, ms->numSimulations_);
117     buf[ms->simulationIndex_] = n;
118     gmx_sumi_sim(ms->numSimulations_, buf, ms);
119     bPos   = TRUE;
120     bEqual = TRUE;
121     for (s = 0; s < ms->numSimulations_; s++)
122     {
123         bPos   = bPos && (buf[s] > 0);
124         bEqual = bEqual && (buf[s] == buf[0]);
125     }
126     if (bPos)
127     {
128         if (bEqual)
129         {
130             nmin = std::min(nmin, buf[0]);
131         }
132         else
133         {
134             /* Find the least common multiple */
135             for (d = 2; d < nmin; d++)
136             {
137                 s = 0;
138                 while (s < ms->numSimulations_ && d % buf[s] == 0)
139                 {
140                     s++;
141                 }
142                 if (s == ms->numSimulations_)
143                 {
144                     /* We found the LCM and it is less than nmin */
145                     nmin = d;
146                     break;
147                 }
148             }
149         }
150     }
151     sfree(buf);
152
153     return nmin;
154 }
155
156 static void calc_ke_part_normal(gmx::ArrayRef<const gmx::RVec> v,
157                                 const t_grpopts*               opts,
158                                 const t_mdatoms*               md,
159                                 gmx_ekindata_t*                ekind,
160                                 t_nrnb*                        nrnb,
161                                 gmx_bool                       bEkinAveVel)
162 {
163     int                         g;
164     gmx::ArrayRef<t_grp_tcstat> tcstat  = ekind->tcstat;
165     gmx::ArrayRef<t_grp_acc>    grpstat = ekind->grpstat;
166
167     /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin.  Leap with AveVel is also
168        an option, but not supported now.
169        bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
170      */
171
172     /* group velocities are calculated in update_ekindata and
173      * accumulated in acumulate_groups.
174      * Now the partial global and groups ekin.
175      */
176     for (g = 0; (g < opts->ngtc); g++)
177     {
178         copy_mat(tcstat[g].ekinh, tcstat[g].ekinh_old);
179         if (bEkinAveVel)
180         {
181             clear_mat(tcstat[g].ekinf);
182             tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
183         }
184         else
185         {
186             clear_mat(tcstat[g].ekinh);
187         }
188     }
189     ekind->dekindl_old = ekind->dekindl;
190     int nthread        = gmx_omp_nthreads_get(emntUpdate);
191
192 #pragma omp parallel for num_threads(nthread) schedule(static)
193     for (int thread = 0; thread < nthread; thread++)
194     {
195         // This OpenMP only loops over arrays and does not call any functions
196         // or memory allocation. It should not be able to throw, so for now
197         // we do not need a try/catch wrapper.
198         int     start_t, end_t, n;
199         int     ga, gt;
200         rvec    v_corrt;
201         real    hm;
202         int     d, m;
203         matrix* ekin_sum;
204         real*   dekindl_sum;
205
206         start_t = ((thread + 0) * md->homenr) / nthread;
207         end_t   = ((thread + 1) * md->homenr) / nthread;
208
209         ekin_sum    = ekind->ekin_work[thread];
210         dekindl_sum = ekind->dekindl_work[thread];
211
212         for (gt = 0; gt < opts->ngtc; gt++)
213         {
214             clear_mat(ekin_sum[gt]);
215         }
216         *dekindl_sum = 0.0;
217
218         ga = 0;
219         gt = 0;
220         for (n = start_t; n < end_t; n++)
221         {
222             if (md->cACC)
223             {
224                 ga = md->cACC[n];
225             }
226             if (md->cTC)
227             {
228                 gt = md->cTC[n];
229             }
230             hm = 0.5 * md->massT[n];
231
232             for (d = 0; (d < DIM); d++)
233             {
234                 v_corrt[d] = v[n][d] - grpstat[ga].u[d];
235             }
236             for (d = 0; (d < DIM); d++)
237             {
238                 for (m = 0; (m < DIM); m++)
239                 {
240                     /* if we're computing a full step velocity, v_corrt[d] has v(t).  Otherwise, v(t+dt/2) */
241                     ekin_sum[gt][m][d] += hm * v_corrt[m] * v_corrt[d];
242                 }
243             }
244             if (md->nMassPerturbed && md->bPerturbed[n])
245             {
246                 *dekindl_sum += 0.5 * (md->massB[n] - md->massA[n]) * iprod(v_corrt, v_corrt);
247             }
248         }
249     }
250
251     ekind->dekindl = 0;
252     for (int thread = 0; thread < nthread; thread++)
253     {
254         for (g = 0; g < opts->ngtc; g++)
255         {
256             if (bEkinAveVel)
257             {
258                 m_add(tcstat[g].ekinf, ekind->ekin_work[thread][g], tcstat[g].ekinf);
259             }
260             else
261             {
262                 m_add(tcstat[g].ekinh, ekind->ekin_work[thread][g], tcstat[g].ekinh);
263             }
264         }
265
266         ekind->dekindl += *ekind->dekindl_work[thread];
267     }
268
269     inc_nrnb(nrnb, eNR_EKIN, md->homenr);
270 }
271
272 static void calc_ke_part_visc(const matrix                   box,
273                               gmx::ArrayRef<const gmx::RVec> x,
274                               gmx::ArrayRef<const gmx::RVec> v,
275                               const t_grpopts*               opts,
276                               const t_mdatoms*               md,
277                               gmx_ekindata_t*                ekind,
278                               t_nrnb*                        nrnb,
279                               gmx_bool                       bEkinAveVel)
280 {
281     int                         start = 0, homenr = md->homenr;
282     int                         g, d, n, m, gt = 0;
283     rvec                        v_corrt;
284     real                        hm;
285     gmx::ArrayRef<t_grp_tcstat> tcstat = ekind->tcstat;
286     t_cos_acc*                  cosacc = &(ekind->cosacc);
287     real                        dekindl;
288     real                        fac, cosz;
289     double                      mvcos;
290
291     for (g = 0; g < opts->ngtc; g++)
292     {
293         copy_mat(ekind->tcstat[g].ekinh, ekind->tcstat[g].ekinh_old);
294         clear_mat(ekind->tcstat[g].ekinh);
295     }
296     ekind->dekindl_old = ekind->dekindl;
297
298     fac     = 2 * M_PI / box[ZZ][ZZ];
299     mvcos   = 0;
300     dekindl = 0;
301     for (n = start; n < start + homenr; n++)
302     {
303         if (md->cTC)
304         {
305             gt = md->cTC[n];
306         }
307         hm = 0.5 * md->massT[n];
308
309         /* Note that the times of x and v differ by half a step */
310         /* MRS -- would have to be changed for VV */
311         cosz = std::cos(fac * x[n][ZZ]);
312         /* Calculate the amplitude of the new velocity profile */
313         mvcos += 2 * cosz * md->massT[n] * v[n][XX];
314
315         copy_rvec(v[n], v_corrt);
316         /* Subtract the profile for the kinetic energy */
317         v_corrt[XX] -= cosz * cosacc->vcos;
318         for (d = 0; (d < DIM); d++)
319         {
320             for (m = 0; (m < DIM); m++)
321             {
322                 /* if we're computing a full step velocity, v_corrt[d] has v(t).  Otherwise, v(t+dt/2) */
323                 if (bEkinAveVel)
324                 {
325                     tcstat[gt].ekinf[m][d] += hm * v_corrt[m] * v_corrt[d];
326                 }
327                 else
328                 {
329                     tcstat[gt].ekinh[m][d] += hm * v_corrt[m] * v_corrt[d];
330                 }
331             }
332         }
333         if (md->nPerturbed && md->bPerturbed[n])
334         {
335             /* The minus sign here might be confusing.
336              * The kinetic contribution from dH/dl doesn't come from
337              * d m(l)/2 v^2 / dl, but rather from d p^2/2m(l) / dl,
338              * where p are the momenta. The difference is only a minus sign.
339              */
340             dekindl -= 0.5 * (md->massB[n] - md->massA[n]) * iprod(v_corrt, v_corrt);
341         }
342     }
343     ekind->dekindl = dekindl;
344     cosacc->mvcos  = mvcos;
345
346     inc_nrnb(nrnb, eNR_EKIN, homenr);
347 }
348
349 static void calc_ke_part(gmx::ArrayRef<const gmx::RVec> x,
350                          gmx::ArrayRef<const gmx::RVec> v,
351                          const matrix                   box,
352                          const t_grpopts*               opts,
353                          const t_mdatoms*               md,
354                          gmx_ekindata_t*                ekind,
355                          t_nrnb*                        nrnb,
356                          gmx_bool                       bEkinAveVel)
357 {
358     if (ekind->cosacc.cos_accel == 0)
359     {
360         calc_ke_part_normal(v, opts, md, ekind, nrnb, bEkinAveVel);
361     }
362     else
363     {
364         calc_ke_part_visc(box, x, v, opts, md, ekind, nrnb, bEkinAveVel);
365     }
366 }
367
368 /* TODO Specialize this routine into init-time and loop-time versions?
369    e.g. bReadEkin is only true when restoring from checkpoint */
370 void compute_globals(gmx_global_stat*               gstat,
371                      t_commrec*                     cr,
372                      const t_inputrec*              ir,
373                      t_forcerec*                    fr,
374                      gmx_ekindata_t*                ekind,
375                      gmx::ArrayRef<const gmx::RVec> x,
376                      gmx::ArrayRef<const gmx::RVec> v,
377                      const matrix                   box,
378                      const t_mdatoms*               mdatoms,
379                      t_nrnb*                        nrnb,
380                      t_vcm*                         vcm,
381                      gmx_wallcycle_t                wcycle,
382                      gmx_enerdata_t*                enerd,
383                      tensor                         force_vir,
384                      tensor                         shake_vir,
385                      tensor                         total_vir,
386                      tensor                         pres,
387                      gmx::Constraints*              constr,
388                      gmx::SimulationSignaller*      signalCoordinator,
389                      const matrix                   lastbox,
390                      int*                           totalNumberOfBondedInteractions,
391                      gmx_bool*                      bSumEkinhOld,
392                      const int                      flags)
393 {
394     gmx_bool bEner, bPres, bTemp;
395     gmx_bool bStopCM, bGStat, bReadEkin, bEkinAveVel, bScaleEkin, bConstrain;
396     gmx_bool bCheckNumberOfBondedInteractions;
397     real     dvdl_ekin;
398
399     /* translate CGLO flags to gmx_booleans */
400     bStopCM                          = ((flags & CGLO_STOPCM) != 0);
401     bGStat                           = ((flags & CGLO_GSTAT) != 0);
402     bReadEkin                        = ((flags & CGLO_READEKIN) != 0);
403     bScaleEkin                       = ((flags & CGLO_SCALEEKIN) != 0);
404     bEner                            = ((flags & CGLO_ENERGY) != 0);
405     bTemp                            = ((flags & CGLO_TEMPERATURE) != 0);
406     bPres                            = ((flags & CGLO_PRESSURE) != 0);
407     bConstrain                       = ((flags & CGLO_CONSTRAINT) != 0);
408     bCheckNumberOfBondedInteractions = ((flags & CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS) != 0);
409
410     /* we calculate a full state kinetic energy either with full-step velocity verlet
411        or half step where we need the pressure */
412
413     bEkinAveVel = (ir->eI == eiVV || (ir->eI == eiVVAK && bPres) || bReadEkin);
414
415     /* in initalization, it sums the shake virial in vv, and to
416        sums ekinh_old in leapfrog (or if we are calculating ekinh_old) for other reasons */
417
418     /* ########## Kinetic energy  ############## */
419
420     if (bTemp)
421     {
422         /* Non-equilibrium MD: this is parallellized, but only does communication
423          * when there really is NEMD.
424          */
425
426         if (PAR(cr) && (ekind->bNEMD))
427         {
428             accumulate_u(cr, &(ir->opts), ekind);
429         }
430         if (!bReadEkin)
431         {
432             calc_ke_part(x, v, box, &(ir->opts), mdatoms, ekind, nrnb, bEkinAveVel);
433         }
434     }
435
436     /* Calculate center of mass velocity if necessary, also parallellized */
437     if (bStopCM)
438     {
439         calc_vcm_grp(*mdatoms, x, v, vcm);
440     }
441
442     if (bTemp || bStopCM || bPres || bEner || bConstrain || bCheckNumberOfBondedInteractions)
443     {
444         if (!bGStat)
445         {
446             /* We will not sum ekinh_old,
447              * so signal that we still have to do it.
448              */
449             *bSumEkinhOld = TRUE;
450         }
451         else
452         {
453             gmx::ArrayRef<real> signalBuffer = signalCoordinator->getCommunicationBuffer();
454             if (PAR(cr))
455             {
456                 wallcycle_start(wcycle, ewcMoveE);
457                 global_stat(gstat, cr, enerd, force_vir, shake_vir, ir, ekind, constr,
458                             bStopCM ? vcm : nullptr, signalBuffer.size(), signalBuffer.data(),
459                             totalNumberOfBondedInteractions, *bSumEkinhOld, flags);
460                 wallcycle_stop(wcycle, ewcMoveE);
461             }
462             signalCoordinator->finalizeSignals();
463             *bSumEkinhOld = FALSE;
464         }
465     }
466
467     if (bEner)
468     {
469         /* Calculate the amplitude of the cosine velocity profile */
470         ekind->cosacc.vcos = ekind->cosacc.mvcos / mdatoms->tmass;
471     }
472
473     if (bTemp)
474     {
475         /* Sum the kinetic energies of the groups & calc temp */
476         /* compute full step kinetic energies if vv, or if vv-avek and we are computing the pressure with inputrecNptTrotter */
477         /* three maincase:  VV with AveVel (md-vv), vv with AveEkin (md-vv-avek), leap with AveEkin (md).
478            Leap with AveVel is not supported; it's not clear that it will actually work.
479            bEkinAveVel: If TRUE, we simply multiply ekin by ekinscale to get a full step kinetic energy.
480            If FALSE, we average ekinh_old and ekinh*ekinscale_nhc to get an averaged half step kinetic energy.
481          */
482         enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, &dvdl_ekin, bEkinAveVel, bScaleEkin);
483         enerd->dvdl_lin[efptMASS] = static_cast<double>(dvdl_ekin);
484
485         enerd->term[F_EKIN] = trace(ekind->ekin);
486
487         for (auto& dhdl : enerd->dhdlLambda)
488         {
489             dhdl += enerd->dvdl_lin[efptMASS];
490         }
491     }
492
493     /* ########## Now pressure ############## */
494     // TODO: For the VV integrator bConstrain is needed in the conditional. This is confusing, so get rid of this.
495     if (bPres || bConstrain)
496     {
497         m_add(force_vir, shake_vir, total_vir);
498
499         /* Calculate pressure and apply LR correction if PPPM is used.
500          * Use the box from last timestep since we already called update().
501          */
502
503         enerd->term[F_PRES] = calc_pres(fr->pbcType, ir->nwall, lastbox, ekind->ekin, total_vir, pres);
504     }
505 }
506
507 void setCurrentLambdasRerun(int64_t           step,
508                             const t_lambda*   fepvals,
509                             const t_trxframe* rerun_fr,
510                             const double*     lam0,
511                             t_state*          globalState)
512 {
513     GMX_RELEASE_ASSERT(globalState != nullptr,
514                        "setCurrentLambdasGlobalRerun should be called with a valid state object");
515
516     if (rerun_fr->bLambda)
517     {
518         if (fepvals->delta_lambda == 0)
519         {
520             globalState->lambda[efptFEP] = rerun_fr->lambda;
521         }
522         else
523         {
524             /* find out between which two value of lambda we should be */
525             real frac      = step * fepvals->delta_lambda;
526             int  fep_state = static_cast<int>(std::floor(frac * fepvals->n_lambda));
527             /* interpolate between this state and the next */
528             /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
529             frac = frac * fepvals->n_lambda - fep_state;
530             for (int i = 0; i < efptNR; i++)
531             {
532                 globalState->lambda[i] =
533                         lam0[i] + (fepvals->all_lambda[i][fep_state])
534                         + frac * (fepvals->all_lambda[i][fep_state + 1] - fepvals->all_lambda[i][fep_state]);
535             }
536         }
537     }
538     else if (rerun_fr->bFepState)
539     {
540         globalState->fep_state = rerun_fr->fep_state;
541         for (int i = 0; i < efptNR; i++)
542         {
543             globalState->lambda[i] = fepvals->all_lambda[i][globalState->fep_state];
544         }
545     }
546 }
547
548 void setCurrentLambdasLocal(const int64_t       step,
549                             const t_lambda*     fepvals,
550                             const double*       lam0,
551                             gmx::ArrayRef<real> lambda,
552                             const int           currentFEPState)
553 /* find the current lambdas.  If rerunning, we either read in a state, or a lambda value,
554    requiring different logic. */
555 {
556     if (fepvals->delta_lambda != 0)
557     {
558         /* find out between which two value of lambda we should be */
559         real frac = step * fepvals->delta_lambda;
560         if (fepvals->n_lambda > 0)
561         {
562             int fep_state = static_cast<int>(std::floor(frac * fepvals->n_lambda));
563             /* interpolate between this state and the next */
564             /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
565             frac = frac * fepvals->n_lambda - fep_state;
566             for (int i = 0; i < efptNR; i++)
567             {
568                 lambda[i] = lam0[i] + (fepvals->all_lambda[i][fep_state])
569                             + frac * (fepvals->all_lambda[i][fep_state + 1] - fepvals->all_lambda[i][fep_state]);
570             }
571         }
572         else
573         {
574             for (int i = 0; i < efptNR; i++)
575             {
576                 lambda[i] = lam0[i] + frac;
577             }
578         }
579     }
580     else
581     {
582         /* if < 0, fep_state was never defined, and we should not set lambda from the state */
583         if (currentFEPState > -1)
584         {
585             for (int i = 0; i < efptNR; i++)
586             {
587                 lambda[i] = fepvals->all_lambda[i][currentFEPState];
588             }
589         }
590     }
591 }
592
593 static void min_zero(int* n, int i)
594 {
595     if (i > 0 && (*n == 0 || i < *n))
596     {
597         *n = i;
598     }
599 }
600
601 static int lcd4(int i1, int i2, int i3, int i4)
602 {
603     int nst;
604
605     nst = 0;
606     min_zero(&nst, i1);
607     min_zero(&nst, i2);
608     min_zero(&nst, i3);
609     min_zero(&nst, i4);
610     if (nst == 0)
611     {
612         gmx_incons("All 4 inputs for determining nstglobalcomm are <= 0");
613     }
614
615     while (nst > 1
616            && ((i1 > 0 && i1 % nst != 0) || (i2 > 0 && i2 % nst != 0) || (i3 > 0 && i3 % nst != 0)
617                || (i4 > 0 && i4 % nst != 0)))
618     {
619         nst--;
620     }
621
622     return nst;
623 }
624
625 int computeGlobalCommunicationPeriod(const gmx::MDLogger& mdlog, t_inputrec* ir, const t_commrec* cr)
626 {
627     int nstglobalcomm;
628     {
629         // Set up the default behaviour
630         if (!(ir->nstcalcenergy > 0 || ir->nstlist > 0 || ir->etc != etcNO || ir->epc != epcNO))
631         {
632             /* The user didn't choose the period for anything
633                important, so we just make sure we can send signals and
634                write output suitably. */
635             nstglobalcomm = 10;
636             if (ir->nstenergy > 0 && ir->nstenergy < nstglobalcomm)
637             {
638                 nstglobalcomm = ir->nstenergy;
639             }
640         }
641         else
642         {
643             /* The user has made a choice (perhaps implicitly), so we
644              * ensure that we do timely intra-simulation communication
645              * for (possibly) each of the four parts that care.
646              *
647              * TODO Does the Verlet scheme (+ DD) need any
648              * communication at nstlist steps? Is the use of nstlist
649              * here a leftover of the twin-range scheme? Can we remove
650              * nstlist when we remove the group scheme?
651              */
652             nstglobalcomm = lcd4(ir->nstcalcenergy, ir->nstlist, ir->etc != etcNO ? ir->nsttcouple : 0,
653                                  ir->epc != epcNO ? ir->nstpcouple : 0);
654         }
655     }
656
657     // TODO change this behaviour. Instead grompp should print
658     // a (performance) note and mdrun should not change ir.
659     if (ir->comm_mode != ecmNO && ir->nstcomm < nstglobalcomm)
660     {
661         GMX_LOG(mdlog.warning)
662                 .asParagraph()
663                 .appendTextFormatted("WARNING: Changing nstcomm from %d to %d", ir->nstcomm, nstglobalcomm);
664         ir->nstcomm = nstglobalcomm;
665     }
666
667     if (cr->nnodes > 1)
668     {
669         GMX_LOG(mdlog.info)
670                 .appendTextFormatted("Intra-simulation communication will occur every %d steps.\n",
671                                      nstglobalcomm);
672     }
673     return nstglobalcomm;
674 }
675
676 void rerun_parallel_comm(t_commrec* cr, t_trxframe* fr, gmx_bool* bLastStep)
677 {
678     rvec *xp, *vp;
679
680     if (MASTER(cr) && *bLastStep)
681     {
682         fr->natoms = -1;
683     }
684     xp = fr->x;
685     vp = fr->v;
686     gmx_bcast(sizeof(*fr), fr, cr->mpi_comm_mygroup);
687     fr->x = xp;
688     fr->v = vp;
689
690     *bLastStep = (fr->natoms < 0);
691 }
692
693 // TODO Most of this logic seems to belong in the respective modules
694 void set_state_entries(t_state* state, const t_inputrec* ir, bool useModularSimulator)
695 {
696     /* The entries in the state in the tpx file might not correspond
697      * with what is needed, so we correct this here.
698      */
699     state->flags = 0;
700     if (ir->efep != efepNO || ir->bExpanded)
701     {
702         state->flags |= (1 << estLAMBDA);
703         state->flags |= (1 << estFEPSTATE);
704     }
705     state->flags |= (1 << estX);
706     GMX_RELEASE_ASSERT(state->x.size() == state->natoms,
707                        "We should start a run with an initialized state->x");
708     if (EI_DYNAMICS(ir->eI))
709     {
710         state->flags |= (1 << estV);
711     }
712
713     state->nnhpres = 0;
714     if (ir->pbcType != PbcType::No)
715     {
716         state->flags |= (1 << estBOX);
717         if (inputrecPreserveShape(ir))
718         {
719             state->flags |= (1 << estBOX_REL);
720         }
721         if ((ir->epc == epcPARRINELLORAHMAN) || (ir->epc == epcMTTK))
722         {
723             state->flags |= (1 << estBOXV);
724             if (!useModularSimulator)
725             {
726                 state->flags |= (1 << estPRES_PREV);
727             }
728         }
729         if (inputrecNptTrotter(ir) || (inputrecNphTrotter(ir)))
730         {
731             state->nnhpres = 1;
732             state->flags |= (1 << estNHPRES_XI);
733             state->flags |= (1 << estNHPRES_VXI);
734             state->flags |= (1 << estSVIR_PREV);
735             state->flags |= (1 << estFVIR_PREV);
736             state->flags |= (1 << estVETA);
737             state->flags |= (1 << estVOL0);
738         }
739         if (ir->epc == epcBERENDSEN)
740         {
741             state->flags |= (1 << estBAROS_INT);
742         }
743     }
744
745     if (ir->etc == etcNOSEHOOVER)
746     {
747         state->flags |= (1 << estNH_XI);
748         state->flags |= (1 << estNH_VXI);
749     }
750
751     if (ir->etc == etcVRESCALE || ir->etc == etcBERENDSEN)
752     {
753         state->flags |= (1 << estTHERM_INT);
754     }
755
756     init_gtc_state(state, state->ngtc, state->nnhpres,
757                    ir->opts.nhchainlength); /* allocate the space for nose-hoover chains */
758     init_ekinstate(&state->ekinstate, ir);
759
760     if (ir->bExpanded)
761     {
762         snew(state->dfhist, 1);
763         init_df_history(state->dfhist, ir->fepvals->n_lambda);
764     }
765
766     if (ir->pull && ir->pull->bSetPbcRefToPrevStepCOM)
767     {
768         state->flags |= (1 << estPULLCOMPREVSTEP);
769     }
770 }